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 //-------------------------------------------------------------------------
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  if (m_saveSummaryPlot)
53  m_plotter = new TrackingMaterialPlotter(300., 120., 10); // 10x10 points per cm2
54  else
55  m_plotter = nullptr;
56 }
57 
58 //-------------------------------------------------------------------------
60  if (m_plotter)
61  delete m_plotter;
62 }
63 
64 //-------------------------------------------------------------------------
66  std::ofstream parameters(name);
67  std::cout << std::endl;
68  for (unsigned int i = 0; i < m_groups.size(); ++i) {
69  MaterialAccountingGroup& layer = *(m_groups[i]);
70  std::cout << layer.name() << std::endl;
71  std::cout << boost::format("\tnumber of hits: %9d") % layer.tracks() << std::endl;
72  std::cout << boost::format("\tnormalized segment length: %9.1f ± %9.1f cm") % layer.averageLength() %
73  layer.sigmaLength()
74  << std::endl;
75  std::cout << boost::format("\tnormalized radiation lengths: %9.3f ± %9.3f") % layer.averageRadiationLengths() %
76  layer.sigmaRadiationLengths()
77  << std::endl;
78  std::cout << boost::format("\tnormalized energy loss: %6.5fe-03 ± %6.5fe-03 GeV") %
79  layer.averageEnergyLoss() % layer.sigmaEnergyLoss()
80  << std::endl;
81  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") %
82  layer.name() % layer.tracks() % layer.averageLength() % layer.sigmaLength() %
84  layer.sigmaEnergyLoss()
85  << std::endl;
86  }
87  std::cout << std::endl;
88 
89  parameters.close();
90 }
91 
92 //-------------------------------------------------------------------------
94  std::ofstream xml(name);
95  xml << "<?xml version=\"1.0\" encoding=\"utf-8\"?>" << std::endl;
96  xml << "<Groups>" << std::endl;
97  for (unsigned int i = 0; i < m_groups.size(); ++i) {
98  MaterialAccountingGroup& layer = *(m_groups[i]);
99  xml << " <Group name=\"" << layer.name() << "\">\n"
100  << " <Parameter name=\"TrackerRadLength\" value=\"" << layer.averageRadiationLengths() << "\"/>\n"
101  << " <Parameter name=\"TrackerXi\" value=\"" << layer.averageEnergyLoss() << "\"/>\n"
102  << " </Group>\n"
103  << std::endl;
104  }
105  xml << "</Groups>" << std::endl;
106 }
107 
108 //-------------------------------------------------------------------------
110  for (unsigned int i = 0; i < m_groups.size(); ++i) {
111  MaterialAccountingGroup& layer = *(m_groups[i]);
112  layer.savePlots();
113  }
114 }
115 
116 //-------------------------------------------------------------------------
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  using namespace edm;
139  setup.get<IdealGeometryRecord>().get(hDDD);
140 
141  m_groups.reserve(m_groupNames.size());
142  // Initialize m_groups iff it has size equal to zero, so that we are
143  // sure it will never be repopulated with the same entries over and
144  // over again in the eventloop, at each call of the analyze method.
145  if (m_groups.empty()) {
146  for (unsigned int i = 0; i < m_groupNames.size(); ++i)
147  m_groups.push_back(new MaterialAccountingGroup(m_groupNames[i], *hDDD));
148 
149  LogInfo("TrackingMaterialAnalyser") << "TrackingMaterialAnalyser: List of the tracker groups: " << std::endl;
150  for (unsigned int i = 0; i < m_groups.size(); ++i)
151  LogInfo("TrackingMaterialAnalyser") << i << " TrackingMaterialAnalyser:\t" << m_groups[i]->info() << std::endl;
152  }
154  event.getByToken(m_materialToken, h_tracks);
155 
156  for (std::vector<MaterialAccountingTrack>::const_iterator t = h_tracks->begin(), end = h_tracks->end(); t != end;
157  ++t) {
159  split(track);
160  }
161 }
162 
163 //-------------------------------------------------------------------
164 // split a track in segments, each associated to a sensitive detector
165 // in a DetLayer; then, associate each step to one segment, splitting
166 // the steps across the segment boundaries
167 //
168 // Nota Bene: this implementation assumes that the steps stored along
169 // each track are consecutive and adjacent, and that no step can span
170 // across 3 layers, since all steps should split at layer boundaries
171 
173  using namespace edm;
174  // group sensitive detectors by their DetLayer
175  std::vector<int> group(track.detectors().size());
176  for (unsigned int i = 0; i < track.detectors().size(); ++i)
177  group[i] = findLayer(track.detectors()[i]);
178 
179  for (unsigned int i = 0; i < group.size(); ++i)
180  if (group[i] > 0)
181  LogInfo("TrackingMaterialAnalyser") << "For detector i: " << i << " index: " << group[i]
182  << " R-ranges: " << m_groups[group[i] - 1]->getBoundingR().first << ", "
183  << m_groups[group[i] - 1]->getBoundingR().second << group[i]
184  << " Z-ranges: " << m_groups[group[i] - 1]->getBoundingZ().first << ", "
185  << m_groups[group[i] - 1]->getBoundingZ().second << std::endl;
186 
187  unsigned int detectors = track.detectors().size();
188  if (detectors == 0) {
189  // the track doesn't cross any active detector:
190  // keep al material as unassigned
191  if (m_plotter)
192  for (unsigned int i = 1; i < track.steps().size(); ++i)
194  } else {
195  const double TOLERANCE = 0.0001; // 1 um tolerance
196  std::vector<double> limits(detectors + 2);
197 
198  // define the trivial limits
200  limits[0] = track.detectors()[0].m_curvilinearIn - TOLERANCE;
201  else
202  limits[0] = -TOLERANCE;
204  limits[detectors] = track.detectors()[detectors - 1].m_curvilinearOut + TOLERANCE;
205  else
206  limits[detectors] = track.summary().length() + TOLERANCE;
207  limits[detectors + 1] = INFINITY; // this is probably no more needed, but doesn't harm...
208 
209  // pick the algorithm to define the non-trivial limits
210  switch (m_splitMode) {
211  // assign each segment to the the nearest layer
212  // e.g. the material between pixel barrel 3 and TIB 1 will be split among the two
213  case NEAREST_LAYER:
214  for (unsigned int i = 1; i < detectors; ++i) {
215  limits[i] = (track.detectors()[i - 1].m_curvilinearOut + track.detectors()[i].m_curvilinearIn) / 2.;
216  }
217  break;
218 
219  // assign each segment to the the inner layer
220  // e.g. all material between pixel barrel 3 and TIB 1 will go into the pixel barrel
221  case INNER_LAYER:
222  for (unsigned int i = 1; i < detectors; ++i)
223  limits[i] = track.detectors()[i].m_curvilinearIn - TOLERANCE;
224  break;
225 
226  // assign each segment to the the outer layer
227  // e.g. all material between pixel barrel 3 and TIB 1 will go into the TIB
228  case OUTER_LAYER:
229  for (unsigned int i = 1; i < detectors; ++i)
230  limits[i] = track.detectors()[i - 1].m_curvilinearOut + TOLERANCE;
231  break;
232 
233  case UNDEFINED:
234  default:
235  // throw something
236  throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode";
237  }
238 
239  double begin = 0.; // beginning of step, along the track
240  double end = 0.; // end of step, along the track
241  unsigned int i = 1; // step conter
242 
243  // skip the material before the first layer
244  while (end < limits[0]) {
245  const MaterialAccountingStep& step = track.steps()[i++];
246  end = begin + step.length();
247 
248  // do not account material before the first layer
249  if (m_plotter)
251 
252  begin = end;
253  }
254 
255  unsigned int index = 0; // which detector
256  while (i < track.steps().size()) {
257  const MaterialAccountingStep& step = track.steps()[i++];
258 
259  end = begin + step.length();
260 
261  if (begin > limits[detectors]) {
262  // segment after last layer and skipping requested in configuation
263  if (m_plotter)
265  begin = end;
266  continue;
267  }
268 
269  // from here onwards we should be in the accountable region, either completely in a single layer:
270  // limits[index] <= begin < end <= limits[index+1]
271  // or possibly split between 2 layers
272  // limits[index] < begin < limits[index+1] < end < limits[index+2]
273  if (begin < limits[index] or end > limits[index + 2]) {
274  // sanity check
275  std::cerr << "MaterialAccountingTrack::split(): ERROR: internal logic error, expected " << limits[index]
276  << " < " << begin << " < " << limits[index + 1] << std::endl;
277  break;
278  }
279 
280  if (limits[index] <= begin and end <= limits[index + 1]) {
281  // step completely inside current detector range
282  track.detectors()[index].account(step, begin, end);
283  if (m_plotter)
284  m_plotter->plotSegmentInLayer(step, group[index]);
285  } else {
286  // step shared beteewn two detectors, transition at limits[index+1]
287  double fraction = (limits[index + 1] - begin) / (end - begin);
288  assert(fraction < 1.);
289  std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
290 
291  if (m_plotter) {
292  if (index > 0)
293  m_plotter->plotSegmentInLayer(parts.first, group[index]);
294  else
295  // track outside acceptance, keep as unassocated
296  m_plotter->plotSegmentUnassigned(parts.first);
297 
298  if (index + 1 < detectors)
299  m_plotter->plotSegmentInLayer(parts.second, group[index + 1]);
300  else
301  // track outside acceptance, keep as unassocated
302  m_plotter->plotSegmentUnassigned(parts.second);
303  }
304 
305  track.detectors()[index].account(parts.first, begin, limits[index + 1]);
306  ++index; // next layer
307  if (index < detectors)
308  track.detectors()[index].account(parts.second, limits[index + 1], end);
309  }
310  begin = end;
311  }
312  }
313 
314  // add the material from each detector to its layer (if there is one and only one)
315  for (unsigned int i = 0; i < track.detectors().size(); ++i)
316  if (group[i] != 0)
317  m_groups[group[i] - 1]->addDetector(track.detectors()[i]);
318 
319  // end of track: commit internal buffers and reset the m_groups internal state for a new track
320  for (unsigned int i = 0; i < m_groups.size(); ++i)
321  m_groups[i]->endOfTrack();
322 }
323 
324 //-------------------------------------------------------------------------
325 // find the layer index (0: none, 1-3: PixelBarrel,
326 // 4-7: TID, 8-13: TOB, 14-15,28-29: PixelEndcap,
327 // 16-18,30-32: TID, 19-27,33-41: TEC)
329  int index = 0;
330  size_t inside = 0;
331  for (size_t i = 0; i < m_groups.size(); ++i)
332  if (m_groups[i]->inside(detector)) {
333  ++inside;
334  index = i + 1;
335  }
336  if (inside == 0) {
337  index = 0;
338  std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector does not belong to any DetLayer"
339  << std::endl;
340  std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
341  << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
342  << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
343  << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")" << std::endl;
344  }
345  if (inside > 1) {
346  index = 0;
347  std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector belongs to " << inside << " DetLayers"
348  << std::endl;
349  std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
350  << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
351  << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
352  << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")" << std::endl;
353  }
354 
355  return index;
356 }
357 
358 //-------------------------------------------------------------------------
359 // 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
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:2
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)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
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 &)
T get() const
Definition: EventSetup.h:71
void plotSegmentInLayer(const MaterialAccountingStep &step, int layer)
step
Definition: StallMonitor.cc:94
const MaterialAccountingStep & summary() const
Definition: event.py:1