CMS 3D CMS Logo

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::one::EDAnalyzer<> edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 TrackingMaterialAnalyser (const edm::ParameterSet &)
 
 ~TrackingMaterialAnalyser () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Types

enum  SplitMode { NEAREST_LAYER, INNER_LAYER, OUTER_LAYER, UNDEFINED }
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 
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
 
bool m_isHFNose
 
bool m_isHGCal
 
edm::EDGetTokenT< std::vector< MaterialAccountingTrack > > m_materialToken
 
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::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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)
explicit

Definition at line 29 of file TrackingMaterialAnalyser.cc.

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

29  {
31  consumes<std::vector<MaterialAccountingTrack> >(iPSet.getParameter<edm::InputTag>("MaterialAccounting"));
32  m_groupNames = iPSet.getParameter<std::vector<std::string> >("Groups");
33  const std::string& splitmode = iPSet.getParameter<std::string>("SplitMode");
34  if (strcasecmp(splitmode.c_str(), "NearestLayer") == 0) {
36  } else if (strcasecmp(splitmode.c_str(), "InnerLayer") == 0) {
38  } else if (strcasecmp(splitmode.c_str(), "OuterLayer") == 0) {
40  } else {
43  << "Invalid SplitMode \"" << splitmode
44  << "\". Acceptable values are \"NearestLayer\", \"InnerLayer\", \"OuterLayer\".";
45  }
46  m_skipAfterLastDetector = iPSet.getParameter<bool>("SkipAfterLastDetector");
47  m_skipBeforeFirstDetector = iPSet.getParameter<bool>("SkipBeforeFirstDetector");
48  m_saveSummaryPlot = iPSet.getParameter<bool>("SaveSummaryPlot");
49  m_saveDetailedPlots = iPSet.getParameter<bool>("SaveDetailedPlots");
50  m_saveParameters = iPSet.getParameter<bool>("SaveParameters");
51  m_saveXml = iPSet.getParameter<bool>("SaveXML");
52  m_isHGCal = iPSet.getParameter<bool>("isHGCal");
53  m_isHFNose = iPSet.getParameter<bool>("isHFNose");
54  if (m_saveSummaryPlot) {
55  if (m_isHGCal) {
56  m_plotter = new TrackingMaterialPlotter(550., 300., 10);
57  } else if (m_isHFNose) {
58  m_plotter = new TrackingMaterialPlotter(1200., 350., 10);
59  } else {
60  m_plotter = new TrackingMaterialPlotter(300., 120., 10);
61  } // 10x10 points per cm2
62  } else {
63  m_plotter = nullptr;
64  }
65 }
T getParameter(std::string const &) const
edm::EDGetTokenT< std::vector< MaterialAccountingTrack > > m_materialToken
std::vector< std::string > m_groupNames
TrackingMaterialPlotter * m_plotter
TrackingMaterialAnalyser::~TrackingMaterialAnalyser ( void  )
override

Definition at line 68 of file TrackingMaterialAnalyser.cc.

References m_plotter.

68  {
69  if (m_plotter)
70  delete m_plotter;
71 }
TrackingMaterialPlotter * m_plotter

Member Function Documentation

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

Definition at line 145 of file TrackingMaterialAnalyser.cc.

References end, edm::EventSetup::get(), mps_fire::i, m_groupNames, m_groups, m_materialToken, split(), OrderedSet::t, and HLT_2018_cff::track.

145  {
146  using namespace edm;
148  setup.get<IdealGeometryRecord>().get(hDDD);
149 
150  m_groups.reserve(m_groupNames.size());
151  // Initialize m_groups iff it has size equal to zero, so that we are
152  // sure it will never be repopulated with the same entries over and
153  // over again in the eventloop, at each call of the analyze method.
154  if (m_groups.empty()) {
155  for (unsigned int i = 0; i < m_groupNames.size(); ++i)
156  m_groups.push_back(new MaterialAccountingGroup(m_groupNames[i], *hDDD));
157 
158  LogInfo("TrackingMaterialAnalyser") << "TrackingMaterialAnalyser: List of the tracker groups: " << std::endl;
159  for (unsigned int i = 0; i < m_groups.size(); ++i)
160  LogInfo("TrackingMaterialAnalyser") << i << " TrackingMaterialAnalyser:\t" << m_groups[i]->info() << std::endl;
161  }
163  event.getByToken(m_materialToken, h_tracks);
164 
165  for (std::vector<MaterialAccountingTrack>::const_iterator t = h_tracks->begin(), end = h_tracks->end(); t != end;
166  ++t) {
168  split(track);
169  }
170 }
void split(MaterialAccountingTrack &track)
#define end
Definition: vmac.h:39
edm::EDGetTokenT< std::vector< MaterialAccountingTrack > > m_materialToken
std::vector< std::string > m_groupNames
std::vector< MaterialAccountingGroup * > m_groups
HLT enums.
T get() const
Definition: EventSetup.h:73
void TrackingMaterialAnalyser::beginJob ( void  )
inlineoverrideprivatevirtual
void TrackingMaterialAnalyser::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 126 of file TrackingMaterialAnalyser.cc.

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

Referenced by beginJob().

126  {
127  if (m_saveParameters)
128  saveParameters("parameters");
129 
130  if (m_saveXml)
131  saveXml("parameters.xml");
132 
134  saveLayerPlots();
135 
136  if (m_saveSummaryPlot and m_plotter) {
137  m_plotter->normalize();
138  m_plotter->draw();
139  }
140 }
void saveParameters(const char *name)
TrackingMaterialPlotter * m_plotter
int TrackingMaterialAnalyser::findLayer ( const MaterialAccountingDetector detector)
private

Definition at line 339 of file TrackingMaterialAnalyser.cc.

References beam_dqm_sourceclient-live_cfg::cerr, DEFINE_FWK_MODULE, alignBH_cfg::fixed, mps_fire::i, m_groups, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), MaterialAccountingDetector::position(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by beginJob(), and split().

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

Definition at line 118 of file TrackingMaterialAnalyser.cc.

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

Referenced by beginJob(), and endJob().

118  {
119  for (unsigned int i = 0; i < m_groups.size(); ++i) {
120  MaterialAccountingGroup& layer = *(m_groups[i]);
121  layer.savePlots();
122  }
123 }
std::vector< MaterialAccountingGroup * > m_groups
void savePlots(void)
save the plots
void TrackingMaterialAnalyser::saveParameters ( const char *  name)
private

Definition at line 74 of file TrackingMaterialAnalyser.cc.

References MaterialAccountingGroup::averageEnergyLoss(), MaterialAccountingGroup::averageLength(), MaterialAccountingGroup::averageRadiationLengths(), gather_cfg::cout, dqm-mbProfile::format, mps_fire::i, m_groups, MaterialAccountingGroup::name(), MaterialAccountingGroup::sigmaEnergyLoss(), MaterialAccountingGroup::sigmaLength(), MaterialAccountingGroup::sigmaRadiationLengths(), and MaterialAccountingGroup::tracks().

Referenced by beginJob(), and endJob().

74  {
75  std::ofstream parameters(name);
76  std::cout << std::endl;
77  for (unsigned int i = 0; i < m_groups.size(); ++i) {
78  MaterialAccountingGroup& layer = *(m_groups[i]);
79  std::cout << layer.name() << std::endl;
80  std::cout << boost::format("\tnumber of hits: %9d") % layer.tracks() << std::endl;
81  std::cout << boost::format("\tnormalized segment length: %9.1f ± %9.1f cm") % layer.averageLength() %
82  layer.sigmaLength()
83  << std::endl;
84  std::cout << boost::format("\tnormalized radiation lengths: %9.3f ± %9.3f") % layer.averageRadiationLengths() %
85  layer.sigmaRadiationLengths()
86  << std::endl;
87  std::cout << boost::format("\tnormalized energy loss: %6.5fe-03 ± %6.5fe-03 GeV") %
88  layer.averageEnergyLoss() % layer.sigmaEnergyLoss()
89  << std::endl;
90  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") %
91  layer.name() % layer.tracks() % layer.averageLength() % layer.sigmaLength() %
93  layer.sigmaEnergyLoss()
94  << std::endl;
95  }
96  std::cout << std::endl;
97 
98  parameters.close();
99 }
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
void TrackingMaterialAnalyser::saveXml ( const char *  name)
private

Definition at line 102 of file TrackingMaterialAnalyser.cc.

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

Referenced by beginJob(), and endJob().

102  {
103  std::ofstream xml(name);
104  xml << "<?xml version=\"1.0\" encoding=\"utf-8\"?>" << std::endl;
105  xml << "<Groups>" << std::endl;
106  for (unsigned int i = 0; i < m_groups.size(); ++i) {
107  MaterialAccountingGroup& layer = *(m_groups[i]);
108  xml << " <Group name=\"" << layer.name() << "\">\n"
109  << " <Parameter name=\"TrackerRadLength\" value=\"" << layer.averageRadiationLengths() << "\"/>\n"
110  << " <Parameter name=\"TrackerXi\" value=\"" << layer.averageEnergyLoss() << "\"/>\n"
111  << " </Group>\n"
112  << std::endl;
113  }
114  xml << "</Groups>" << std::endl;
115 }
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 181 of file TrackingMaterialAnalyser.cc.

References begin, beam_dqm_sourceclient-live_cfg::cerr, MuonGeometrySanityCheck_cfi::detectors(), MaterialAccountingTrack::detectors(), end, Exception, findLayer(), HLT_2018_cff::fraction, watchdog::group, mps_fire::i, INNER_LAYER, MaterialAccountingStep::length(), TH2PolyOfflineMaps::limits, edm::errors::LogicError, m_groups, m_plotter, m_skipAfterLastDetector, m_skipBeforeFirstDetector, m_splitMode, NEAREST_LAYER, or, OUTER_LAYER, contentValuesFiles::parts, TrackingMaterialPlotter::plotSegmentInLayer(), TrackingMaterialPlotter::plotSegmentUnassigned(), MaterialAccountingStep::split(), MaterialAccountingTrack::steps(), MaterialAccountingTrack::summary(), and UNDEFINED.

Referenced by analyze(), and beginJob().

181  {
182  using namespace edm;
183  // group sensitive detectors by their DetLayer
184  std::vector<int> group(track.detectors().size());
185  for (unsigned int i = 0; i < track.detectors().size(); ++i)
186  group[i] = findLayer(track.detectors()[i]);
187 
188  for (unsigned int i = 0; i < group.size(); ++i)
189  if (group[i] > 0)
190  LogInfo("TrackingMaterialAnalyser") << "For detector i: " << i << " index: " << group[i]
191  << " R-ranges: " << m_groups[group[i] - 1]->getBoundingR().first << ", "
192  << m_groups[group[i] - 1]->getBoundingR().second << group[i]
193  << " Z-ranges: " << m_groups[group[i] - 1]->getBoundingZ().first << ", "
194  << m_groups[group[i] - 1]->getBoundingZ().second << std::endl;
195 
196  unsigned int detectors = track.detectors().size();
197  if (detectors == 0) {
198  // the track doesn't cross any active detector:
199  // keep al material as unassigned
200  if (m_plotter)
201  for (unsigned int i = 1; i < track.steps().size(); ++i)
203  } else {
204  const double TOLERANCE = 0.0001; // 1 um tolerance
205  std::vector<double> limits(detectors + 2);
206 
207  // define the trivial limits
209  limits[0] = track.detectors()[0].m_curvilinearIn - TOLERANCE;
210  else
211  limits[0] = -TOLERANCE;
213  limits[detectors] = track.detectors()[detectors - 1].m_curvilinearOut + TOLERANCE;
214  else
215  limits[detectors] = track.summary().length() + TOLERANCE;
216  limits[detectors + 1] = INFINITY; // this is probably no more needed, but doesn't harm...
217 
218  // pick the algorithm to define the non-trivial limits
219  switch (m_splitMode) {
220  // assign each segment to the the nearest layer
221  // e.g. the material between pixel barrel 3 and TIB 1 will be split among the two
222  case NEAREST_LAYER:
223  for (unsigned int i = 1; i < detectors; ++i) {
224  limits[i] = (track.detectors()[i - 1].m_curvilinearOut + track.detectors()[i].m_curvilinearIn) / 2.;
225  }
226  break;
227 
228  // assign each segment to the the inner layer
229  // e.g. all material between pixel barrel 3 and TIB 1 will go into the pixel barrel
230  case INNER_LAYER:
231  for (unsigned int i = 1; i < detectors; ++i)
232  limits[i] = track.detectors()[i].m_curvilinearIn - TOLERANCE;
233  break;
234 
235  // assign each segment to the the outer layer
236  // e.g. all material between pixel barrel 3 and TIB 1 will go into the TIB
237  case OUTER_LAYER:
238  for (unsigned int i = 1; i < detectors; ++i)
239  limits[i] = track.detectors()[i - 1].m_curvilinearOut + TOLERANCE;
240  break;
241 
242  case UNDEFINED:
243  [[fallthrough]];
244 
245  default:
246  // throw something
247  throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode";
248  }
249 
250  double begin = 0.; // beginning of step, along the track
251  double end = 0.; // end of step, along the track
252  unsigned int i = 1; // step conter
253 
254  // skip the material before the first layer
255  while (end < limits[0]) {
256  const MaterialAccountingStep& step = track.steps()[i++];
257  end = begin + step.length();
258 
259  // do not account material before the first layer
260  if (m_plotter)
262 
263  begin = end;
264  }
265 
266  unsigned int index = 0; // which detector
267  while (i < track.steps().size()) {
268  const MaterialAccountingStep& step = track.steps()[i++];
269 
270  end = begin + step.length();
271 
272  if (begin > limits[detectors]) {
273  // segment after last layer and skipping requested in configuation
274  if (m_plotter)
276  begin = end;
277  continue;
278  }
279 
280  // from here onwards we should be in the accountable region, either completely in a single layer:
281  // limits[index] <= begin < end <= limits[index+1]
282  // or possibly split between 2 layers
283  // limits[index] < begin < limits[index+1] < end < limits[index+2]
284  if (begin < limits[index] or end > limits[index + 2]) {
285  // sanity check
286  std::cerr << "MaterialAccountingTrack::split(): ERROR: internal logic error, expected " << limits[index]
287  << " < " << begin << " < " << limits[index + 1] << std::endl;
288  break;
289  }
290 
291  if (limits[index] <= begin and end <= limits[index + 1]) {
292  // step completely inside current detector range
293  track.detectors()[index].account(step, begin, end);
294  if (m_plotter)
295  m_plotter->plotSegmentInLayer(step, group[index]);
296  } else {
297  // step shared beteewn two detectors, transition at limits[index+1]
298  double fraction = (limits[index + 1] - begin) / (end - begin);
299  assert(fraction < 1.);
300  std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
301 
302  if (m_plotter) {
303  if (index > 0)
304  m_plotter->plotSegmentInLayer(parts.first, group[index]);
305  else
306  // track outside acceptance, keep as unassocated
307  m_plotter->plotSegmentUnassigned(parts.first);
308 
309  if (index + 1 < detectors)
310  m_plotter->plotSegmentInLayer(parts.second, group[index + 1]);
311  else
312  // track outside acceptance, keep as unassocated
313  m_plotter->plotSegmentUnassigned(parts.second);
314  }
315 
316  track.detectors()[index].account(parts.first, begin, limits[index + 1]);
317  ++index; // next layer
318  if (index < detectors)
319  track.detectors()[index].account(parts.second, limits[index + 1], end);
320  }
321  begin = end;
322  }
323  }
324 
325  // add the material from each detector to its layer (if there is one and only one)
326  for (unsigned int i = 0; i < track.detectors().size(); ++i)
327  if (group[i] != 0)
328  m_groups[group[i] - 1]->addDetector(track.detectors()[i]);
329 
330  // end of track: commit internal buffers and reset the m_groups internal state for a new track
331  for (unsigned int i = 0; i < m_groups.size(); ++i)
332  m_groups[i]->endOfTrack();
333 }
void plotSegmentUnassigned(const MaterialAccountingStep &step)
const std::vector< MaterialAccountingDetector > & detectors() const
const std::vector< MaterialAccountingStep > & steps() const
double length(void) const
int findLayer(const MaterialAccountingDetector &detector)
def detectors(dt=True, csc=True, me42=False, chambers=True, superlayers=False, layers=False)
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::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
#define end
Definition: vmac.h:39
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:32
HLT enums.
void plotSegmentInLayer(const MaterialAccountingStep &step, int layer)
step
Definition: StallMonitor.cc:94
const MaterialAccountingStep & summary() const

Member Data Documentation

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

Definition at line 46 of file TrackingMaterialAnalyser.h.

Referenced by analyze(), and TrackingMaterialAnalyser().

std::vector<MaterialAccountingGroup *> TrackingMaterialAnalyser::m_groups
private
bool TrackingMaterialAnalyser::m_isHFNose
private

Definition at line 44 of file TrackingMaterialAnalyser.h.

Referenced by TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_isHGCal
private

Definition at line 43 of file TrackingMaterialAnalyser.h.

Referenced by TrackingMaterialAnalyser().

edm::EDGetTokenT<std::vector<MaterialAccountingTrack> > TrackingMaterialAnalyser::m_materialToken
private

Definition at line 35 of file TrackingMaterialAnalyser.h.

Referenced by analyze(), and TrackingMaterialAnalyser().

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

Definition at line 40 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_saveParameters
private

Definition at line 41 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_saveSummaryPlot
private

Definition at line 39 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_saveXml
private

Definition at line 42 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_skipAfterLastDetector
private

Definition at line 37 of file TrackingMaterialAnalyser.h.

Referenced by split(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_skipBeforeFirstDetector
private

Definition at line 38 of file TrackingMaterialAnalyser.h.

Referenced by split(), and TrackingMaterialAnalyser().

SplitMode TrackingMaterialAnalyser::m_splitMode
private

Definition at line 36 of file TrackingMaterialAnalyser.h.

Referenced by split(), and TrackingMaterialAnalyser().