CMS 3D CMS Logo

TrackingMaterialAnalyser.cc
Go to the documentation of this file.
1 #include <iostream> // FIXME: switch to MessagLogger & friends
2 #include <fstream>
3 #include <sstream>
4 #include <vector>
5 #include <string>
6 #include <stdexcept>
7 #include <cstring>
8 #include <cstdlib>
9 #include <boost/format.hpp>
10 
16 
21 
27 
28 //-------------------------------------------------------------------------
30 {
31  m_materialToken = consumes<std::vector<MaterialAccountingTrack> >(
32  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 = nullptr;
55 }
56 
57 //-------------------------------------------------------------------------
59 {
60  if (m_plotter)
61  delete m_plotter;
62 }
63 
64 //-------------------------------------------------------------------------
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: %6.5fe-03 ± %6.5fe-03 GeV") % 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 }
88 
89 //-------------------------------------------------------------------------
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 }
105 
106 //-------------------------------------------------------------------------
108 {
109  for (unsigned int i = 0; i < m_groups.size(); ++i) {
110  MaterialAccountingGroup & layer = *(m_groups[i]);
111  layer.savePlots();
112  }
113 }
114 
115 //-------------------------------------------------------------------------
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 }
132 
133 //-------------------------------------------------------------------------
134 
135 //-------------------------------------------------------------------------
137 {
138  using namespace edm;
140  setup.get<IdealGeometryRecord>().get( hDDD );
141 
142  m_groups.reserve( m_groupNames.size() );
143  // Initialize m_groups iff it has size equal to zero, so that we are
144  // sure it will never be repopulated with the same entries over and
145  // over again in the eventloop, at each call of the analyze method.
146  if (m_groups.empty()) {
147  for (unsigned int i = 0; i < m_groupNames.size(); ++i)
148  m_groups.push_back( new MaterialAccountingGroup( m_groupNames[i], * hDDD) );
149 
150  LogInfo("TrackingMaterialAnalyser")
151  << "TrackingMaterialAnalyser: List of the tracker groups: " << std::endl;
152  for (unsigned int i = 0; i < m_groups.size(); ++i)
153  LogInfo("TrackingMaterialAnalyser")
154  << i << " TrackingMaterialAnalyser:\t" << m_groups[i]->info() << std::endl;
155  }
157  event.getByToken(m_materialToken, h_tracks);
158 
159  for (std::vector<MaterialAccountingTrack>::const_iterator t = h_tracks->begin(), end = h_tracks->end(); t != end; ++t) {
161  split( track );
162  }
163 }
164 
165 //-------------------------------------------------------------------
166 // split a track in segments, each associated to a sensitive detector
167 // in a DetLayer; then, associate each step to one segment, splitting
168 // the steps across the segment boundaries
169 //
170 // Nota Bene: this implementation assumes that the steps stored along
171 // each track are consecutive and adjacent, and that no step can span
172 // across 3 layers, since all steps should split at layer boundaries
173 
175 {
176  using namespace edm;
177  // group sensitive detectors by their DetLayer
178  std::vector<int> group( track.detectors().size() );
179  for (unsigned int i = 0; i < track.detectors().size(); ++i)
180  group[i] = findLayer( track.detectors()[i] );
181 
182  for (unsigned int i = 0; i < group.size(); ++i)
183  if (group[i] > 0)
184  LogInfo("TrackingMaterialAnalyser") << "For detector i: " << i << " index: "
185  << group[i] << " R-ranges: "
186  << m_groups[group[i]-1]->getBoundingR().first
187  << ", " << m_groups[group[i]-1]->getBoundingR().second
188  << group[i] << " Z-ranges: "
189  << m_groups[group[i]-1]->getBoundingZ().first
190  << ", " << m_groups[group[i]-1]->getBoundingZ().second
191  << std::endl;
192 
193  unsigned int detectors = track.detectors().size();
194  if (detectors == 0) {
195  // the track doesn't cross any active detector:
196  // keep al material as unassigned
197  if (m_plotter)
198  for (unsigned int i = 1; i < track.steps().size(); ++i)
199  m_plotter->plotSegmentUnassigned( track.steps()[i] );
200  } else {
201  const double TOLERANCE = 0.0001; // 1 um tolerance
202  std::vector<double> limits(detectors + 2);
203 
204  // define the trivial limits
206  limits[0] = track.detectors()[0].m_curvilinearIn - TOLERANCE;
207  else
208  limits[0] = - TOLERANCE;
210  limits[detectors] = track.detectors()[detectors-1].m_curvilinearOut + TOLERANCE;
211  else
212  limits[detectors] = track.summary().length() + TOLERANCE;
213  limits[detectors+1] = INFINITY; // this is probably no more needed, but doesn't harm...
214 
215  // pick the algorithm to define the non-trivial limits
216  switch (m_splitMode) {
217  // assign each segment to the the nearest layer
218  // e.g. the material between pixel barrel 3 and TIB 1 will be split among the two
219  case NEAREST_LAYER:
220  for (unsigned int i = 1; i < detectors; ++i) {
221  limits[i] = (track.detectors()[i-1].m_curvilinearOut + track.detectors()[i].m_curvilinearIn) / 2.;
222  }
223  break;
224 
225  // assign each segment to the the inner layer
226  // e.g. all material between pixel barrel 3 and TIB 1 will go into the pixel barrel
227  case INNER_LAYER:
228  for (unsigned int i = 1; i < detectors; ++i)
229  limits[i] = track.detectors()[i].m_curvilinearIn - TOLERANCE;
230  break;
231 
232  // assign each segment to the the outer layer
233  // e.g. all material between pixel barrel 3 and TIB 1 will go into the TIB
234  case OUTER_LAYER:
235  for (unsigned int i = 1; i < detectors; ++i)
236  limits[i] = track.detectors()[i-1].m_curvilinearOut + TOLERANCE;
237  break;
238 
239  case UNDEFINED:
240  default:
241  // throw something
242  throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode";
243  }
244 
245 
246  double begin = 0.; // beginning of step, along the track
247  double end = 0.; // end of step, along the track
248  unsigned int i = 1; // step conter
249 
250  // skip the material before the first layer
251  while (end < limits[0]) {
252  const MaterialAccountingStep & step = track.steps()[i++];
253  end = begin + step.length();
254 
255  // do not account material before the first layer
256  if (m_plotter)
258 
259  begin = end;
260  }
261 
262  unsigned int index = 0; // which detector
263  while (i < track.steps().size()) {
264  const MaterialAccountingStep & step = track.steps()[i++];
265 
266  end = begin + step.length();
267 
268  if (begin > limits[detectors]) {
269  // segment after last layer and skipping requested in configuation
270  if (m_plotter)
272  begin = end;
273  continue;
274  }
275 
276  // from here onwards we should be in the accountable region, either completely in a single layer:
277  // limits[index] <= begin < end <= limits[index+1]
278  // or possibly split between 2 layers
279  // limits[index] < begin < limits[index+1] < end < limits[index+2]
280  if (begin < limits[index] or end > limits[index+2]) {
281  // sanity check
282  std::cerr << "MaterialAccountingTrack::split(): ERROR: internal logic error, expected " << limits[index] << " < " << begin << " < " << limits[index+1] << std::endl;
283  break;
284  }
285 
286  if (limits[index] <= begin and end <= limits[index+1]) {
287  // step completely inside current detector range
288  track.detectors()[index].account( step, begin, end );
289  if (m_plotter)
290  m_plotter->plotSegmentInLayer( step, group[index] );
291  } else {
292  // step shared beteewn two detectors, transition at limits[index+1]
293  double fraction = (limits[index+1] - begin) / (end - begin);
294  assert(fraction < 1.);
295  std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
296 
297  if (m_plotter) {
298  if (index > 0)
299  m_plotter->plotSegmentInLayer( parts.first, group[index] );
300  else
301  // track outside acceptance, keep as unassocated
302  m_plotter->plotSegmentUnassigned( parts.first );
303 
304  if (index+1 < detectors)
305  m_plotter->plotSegmentInLayer( parts.second, group[index+1] );
306  else
307  // track outside acceptance, keep as unassocated
308  m_plotter->plotSegmentUnassigned( parts.second );
309  }
310 
311  track.detectors()[index].account( parts.first, begin, limits[index+1] );
312  ++index; // next layer
313  if (index < detectors)
314  track.detectors()[index].account( parts.second, limits[index+1], end );
315  }
316  begin = end;
317  }
318 
319  }
320 
321  // add the material from each detector to its layer (if there is one and only one)
322  for (unsigned int i = 0; i < track.detectors().size(); ++i)
323  if (group[i] != 0)
324  m_groups[group[i]-1]->addDetector( track.detectors()[i] );
325 
326  // end of track: commit internal buffers and reset the m_groups internal state for a new track
327  for (unsigned int i = 0; i < m_groups.size(); ++i)
328  m_groups[i]->endOfTrack();
329 }
330 
331 //-------------------------------------------------------------------------
332 // find the layer index (0: none, 1-3: PixelBarrel,
333 // 4-7: TID, 8-13: TOB, 14-15,28-29: PixelEndcap,
334 // 16-18,30-32: TID, 19-27,33-41: TEC)
336 {
337  int index = 0;
338  size_t inside = 0;
339  for (size_t i = 0; i < m_groups.size(); ++i)
340  if (m_groups[i]->inside(detector)) {
341  ++inside;
342  index = i+1;
343  }
344  if (inside == 0) {
345  index = 0;
346  std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector does not belong to any DetLayer" << std::endl;
347  std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
348  << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
349  << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
350  << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")"
351  << std::endl;
352  }
353  if (inside > 1) {
354  index = 0;
355  std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector belongs to " << inside << " DetLayers" << std::endl;
356  std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
357  << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
358  << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
359  << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")"
360  << std::endl;
361  }
362 
363  return index;
364 }
365 
366 //-------------------------------------------------------------------------
367 // define as a plugin
void plotSegmentUnassigned(const MaterialAccountingStep &step)
T getParameter(std::string const &) const
const std::vector< MaterialAccountingDetector > & detectors() const
const std::vector< MaterialAccountingStep > & steps() const
T perp() const
Definition: PV3DBase.h:72
double averageRadiationLengths(void) const
return the average normalized number of radiation lengths
double sigmaLength(void) const
return the sigma of the normalized layer thickness
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double sigmaRadiationLengths(void) const
return the sigma of the normalized number of radiation lengths
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
double length(void) const
double sigmaEnergyLoss(void) const
return the sigma of the normalized energy loss density factor for Bethe-Bloch
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
void split(MaterialAccountingTrack &track)
unsigned int tracks(void) const
return the number of tracks that hit this layer
double averageLength(void) const
return the average normalized layer thickness
void saveXml(const char *name)
void saveParameters(const char *name)
int findLayer(const MaterialAccountingDetector &detector)
def detectors(dt=True, csc=True, me42=False, chambers=True, superlayers=False, layers=False)
double averageEnergyLoss(void) const
return the average normalized energy loss density factor for Bethe-Bloch
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
T z() const
Definition: PV3DBase.h:64
#define end
Definition: vmac.h:39
format
Some error handling for the usage.
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< std::vector< MaterialAccountingTrack > > m_materialToken
const std::string & name(void) const
get the layer name
std::vector< std::string > m_groupNames
TrackingMaterialPlotter * m_plotter
std::vector< MaterialAccountingGroup * > m_groups
const T & get() const
Definition: EventSetup.h:59
void savePlots(void)
save the plots
std::pair< MaterialAccountingStep, MaterialAccountingStep > split(double fraction) const
split the step (0..1) in (0..f) + (f..1) using linear interpolation
const GlobalPoint & position() const
#define begin
Definition: vmac.h:32
HLT enums.
TrackingMaterialAnalyser(const edm::ParameterSet &)
void plotSegmentInLayer(const MaterialAccountingStep &step, int layer)
step
const MaterialAccountingStep & summary() const
Definition: event.py:1