CMS 3D CMS Logo

TrackingRecoMaterialAnalyser.cc
Go to the documentation of this file.
9 
15 
25 
26 #include <cassert>
27 #include <unordered_map>
28 #include <string>
29 
30 // Values are not ordered randomly, but the order is taken from
31 // http://cmslxr.fnal.gov/dxr/CMSSW/source/Geometry/CommonDetUnit/interface/GeomDetEnumerators.h#15
32 static const std::vector<std::string> sDETS{ "", "PXB", "PXF", "TIB", "TID", "TOB", "TEC" };
33 
35  public:
38  edm::Run const&,
39  edm::EventSetup const&) override;
40  void analyze(const edm::Event &, const edm::EventSetup &) override ;
42  private:
43  inline bool isDoubleSided(SiStripDetId strip_id) {
44  return ((strip_id.subDetector() != SiStripDetId::UNKNOWN) && strip_id.glued());
45  }
50  bool usePV_;
52  std::unordered_map<std::string, MonitorElement *> histosOriEta_;
53  std::unordered_map<std::string, MonitorElement *> histosEta_;
67 };
68 
69 //-------------------------------------------------------------------------
71  refitter_(iPSet),
72  tracksToken_(consumes<reco::TrackCollection>(iPSet.getParameter<edm::InputTag>("tracks"))),
73  beamspotToken_(consumes<reco::BeamSpot>(iPSet.getParameter<edm::InputTag>("beamspot"))),
74  verticesToken_(mayConsume<reco::VertexCollection>(iPSet.getParameter<edm::InputTag>("vertices"))),
75  usePV_(iPSet.getParameter<bool>("usePV")),
76  folder_(iPSet.getParameter<std::string>("folder")),
90 {
91 }
92 
93 //-------------------------------------------------------------------------
95 {
96 }
97 
98 //-------------------------------------------------------------------------
100  edm::Run const&,
101  edm::EventSetup const &setup) {
102  using namespace std;
103  edm::ESHandle<TrackerGeometry> trackerGeometry;
104  setup.get<TrackerDigiGeometryRecord>().get(trackerGeometry);
105 
106  ibook.setCurrentFolder(folder_);
107 
108  // Histogram to store the radiation length map, in the R-Z plane,
109  // gathering numbers directly from the trackerRecoMaterial.xml
110  // file. The numbers are not corrected for the track angle.
111  histo_RZ_Ori_ = ibook.bookProfile2D("OriRadLen", "Original_RadLen",
112  600, -300., 300, 120, 0., 120., 0., 1.);
113 
114  // Histogram to store the radiation length map, as before, but
115  // correcting the numbers with the track angle. This represents the
116  // real material seen by the track.
117  histo_RZ_ = ibook.bookProfile2D("RadLen", "RadLen",
118  600, -300., 300, 120, 0., 120., 0., 1.);
119 
120  // Histogram to show the deltaP Out-In in the eta-phi plane.
121  deltaP_in_out_vs_eta_vs_phi_2d_ = ibook.bookProfile2D("DeltaP_in_out_vs_eta_vs_phi_2d",
122  "DeltaP_in_out_vs_eta_vs_phi_2d",
123  100, -3.0, 3.0,
124  100, -3.15, 3.15,
125  0., 100.);
126 
127  // Histogram to show the deltaP Out-In vs eta.
128  deltaP_in_out_vs_eta_2d_ = ibook.book2D("DeltaP_in_out_vs_eta_2d", "DeltaP_in_out_vs_eta_2d",
129  100, -3.0, 3.0, 100, 0., 1);
130 
131  // Histogram to show the deltaP Out-In vs Z. The Z coordinate is the
132  // one computed at the outermost hit.
133  deltaP_in_out_vs_z_2d_ = ibook.book2D("DeltaP_in_out_vs_z_2d", "DeltaP_in_out_vs_z_2d",
134  600, -300, 300, 200., -1, 1.);
135 
136  // Histogram to show the deltaP Out-In vs eta. The eta is the one
137  // computed at the innermost hit.
138  deltaP_in_out_vs_eta_ = ibook.bookProfile("DeltaP_in_out_vs_eta", "DeltaP_in_out_vs_eta",
139  100, -3.0, 3.0, -100., 100.);
140  deltaP_in_out_vs_z_ = ibook.bookProfile("DeltaP_in_out_vs_z", "DeltaP_in_out_vs_z",
141  600, -300, 300, -100., 100.);
142 
143  // Histogram to show the delta_Pt Out-In vs eta. The eta is the one
144  // computed at the innermost hit.
145  deltaPt_in_out_vs_eta_ = ibook.bookProfile("DeltaPt_in_out_vs_eta", "DeltaPt_in_out_vs_eta",
146  100, -3.0, 3.0, -100., 100.);
147 
148  // Histogram to show the delta_Pt Out-In vs Z. The Z is the one
149  // computed at the outermost hit.
150  deltaPt_in_out_vs_z_ = ibook.bookProfile("DeltaPt_in_out_vs_z", "DeltaPt_in_out_vs_z",
151  600, -300, 300, -100., 100);
152 
153  // Histogram to show the delta_Pl Out-In vs eta. The eta is the one
154  // computed at the innermost hit.
155  deltaPl_in_out_vs_eta_ = ibook.bookProfile("DeltaPz_in_out_vs_eta", "DeltaPz_in_out_vs_eta",
156  100, -3.0, 3.0, -100., 100.);
157 
158  // Histogram to show the delta_Pl Out-In vs Z. The Z is the one
159  // computed at the outermost hit.
160  deltaPl_in_out_vs_z_ = ibook.bookProfile("DeltaPz_in_out_vs_z", "DeltaPz_in_out_vs_z",
161  600, -300, 300, -100., 100.);
162 
163  // Histogram to show the delta_Pt Out-In in the Z-R plane. Z and R
164  // are related to the outermost hit.
165  deltaPt_in_out_2d_ = ibook.bookProfile2D("DeltaPt 2D", "DeltaPt 2D",
166  600, -300., 300, 120, 0., 120., -100., 100.);
167 
168  // Histogram to show the distribution of p vs eta for all tracks.
169  P_vs_eta_2d_ = ibook.book2D("P_vs_eta_2d", "P_vs_eta_2d",
170  100, -3.0, 3.0, 100., 0., 5.);
171  char title[50];
172  char key[20];
173  for (unsigned int det = 1; det < sDETS.size(); ++det ) {
174  for (unsigned int sub_det = 1;
175  sub_det <= trackerGeometry->numberOfLayers(det); ++sub_det) {
176  memset(title, 0, sizeof(title));
177  snprintf(title, sizeof(title), "Original_RadLen_vs_Eta_%s%d", sDETS[det].data(), sub_det);
178  snprintf(key, sizeof(key), "%s%d", sDETS[det].data(), sub_det);
179  histosOriEta_.insert(make_pair<string, MonitorElement*>(key,
180  ibook.bookProfile(title, title, 250, -5.0, 5.0, 0., 1.)));
181  snprintf(title, sizeof(title), "RadLen_vs_Eta_%s%d", sDETS[det].data(), sub_det);
182  histosEta_.insert(make_pair<string, MonitorElement*>(key,
183  ibook.bookProfile(title, title, 250, -5.0, 5.0, 0., 1.)));
184  }
185  }
186 }
187 
188 //-------------------------------------------------------------------------
190  const edm::EventSetup& setup)
191 {
192  using namespace edm;
193  using namespace reco;
194  using namespace std;
195 
196  refitter_.setServices(setup);
197 
200  ESHandle<TrackerTopology> trk_topology;
201 
202  // Get the TrackerTopology
203  setup.get<TrackerTopologyRcd>().get(trk_topology);
204 
205  // Get Tracks
206  event.getByToken(tracksToken_, tracks);
207  if (!tracks.isValid() || tracks->empty()) {
208  LogInfo("TrackingRecoMaterialAnalyser") << "Invalid or empty track collection" << endl;
209  return;
210  }
211 
212  // Track Selector
213  auto selector = [](const Track &track, const reco::Vertex::Point &pv) -> bool {
214  return (track.quality(track.qualityByName("highPurity"))
215  && track.dxy(pv) < 0.01
216  && track.hitPattern().numberOfLostTrackerHits(HitPattern::MISSING_OUTER_HITS) == 0);
217  };
218 
219  // Get BeamSpot
221  event.getByToken(beamspotToken_, beamSpot);
222  // Bail out if missing
223  if (!beamSpot.isValid())
224  return;
225 
226  reco::Vertex::Point pv(beamSpot->position());
227  if (usePV_) {
228  event.getByToken(verticesToken_, vertices);
229  if (vertices.isValid() && !vertices->empty()) {
230  // Since we need to use eta and Z information from the tracks, in
231  // order not to have the reco material distribution washed out due
232  // to geometrical effects, we need to confine the PV to some small
233  // region.
234  const Vertex &v = (*vertices)[0];
235  if (!v.isFake() && v.ndof() > 4
236  && std::fabs(v.z()) < 24
237  && std::fabs(v.position().rho()) < 2)
238  pv = v.position();
239  }
240  }
241 
242  // Main idea:
243  // * select first good tracks in input, according to reasonable criteria
244  // * refit the tracks so that we have access to the TrajectoryMeasurements
245  // that internally have all the information about the TSOS on all
246  // crossed layers. We need the refit track and not the original one so
247  // that we are able to correctly compute the path travelled by the track
248  // within the detector, using its updated TSOS. The material description
249  // can in principle be derived also directly from the rechits, via the
250  // det()[(->)GeomDet *]->mediumProperties chain, but that would simply give the
251  // face values, not the "real" ones used while propagating the track.
252  // * Loop on all measurements, extract the information about the TSOS,
253  // its surface and its mediumProperties
254  // * Make plots for the untouched material properties, but also for the
255  // ones corrected by the track direction, since the material properties,
256  // according to the documentation, should refer to normal incidence of
257  // the track, which is seldom the case, according to the current direction
258  TrajectoryStateOnSurface current_tsos;
259  DetId current_det;
260  for (auto const track : *tracks) {
261  if (!selector(track, pv))
262  continue;
263  auto const& inner = track.innerMomentum();
264  auto const& outer = track.outerMomentum();
265  deltaP_in_out_vs_eta_->Fill(inner.eta(), inner.R() - outer.R());
266  deltaP_in_out_vs_z_->Fill(track.outerZ(), inner.R() - outer.R());
267  deltaP_in_out_vs_eta_2d_->Fill(inner.eta(), inner.R() - outer.R());
268  deltaP_in_out_vs_eta_vs_phi_2d_->Fill(inner.eta(), inner.phi(), inner.R() - outer.R());
269  deltaP_in_out_vs_z_2d_->Fill(track.outerZ(), inner.R() - outer.R());
270  deltaPt_in_out_vs_eta_->Fill(inner.eta(), inner.rho() - outer.rho());
271  deltaPt_in_out_vs_z_->Fill(track.outerZ(), inner.rho() - outer.rho());
272  deltaPl_in_out_vs_eta_->Fill(inner.eta(), inner.z() - outer.z());
273  deltaPl_in_out_vs_z_->Fill(track.outerZ(), inner.z() - outer.z());
274  deltaPt_in_out_2d_->Fill(track.outerZ(), track.outerPosition().rho(), inner.rho() - outer.rho());
275  P_vs_eta_2d_->Fill(track.eta(), track.p());
276  vector<Trajectory> traj = refitter_.transform(track);
277  if (traj.size() > 1 || traj.empty())
278  continue;
279  for (auto const &tm : traj.front().measurements()) {
280  if (tm.recHit().get() &&
281  (tm.recHitR().type() == TrackingRecHit::valid ||
282  tm.recHitR().type() == TrackingRecHit::missing)) {
283  current_tsos = tm.updatedState().isValid() ? tm.updatedState() : tm.forwardPredictedState();
284  auto const & localP = current_tsos.localMomentum();
285  current_det = tm.recHit()->geographicalId();
286  const Surface& surface = current_tsos.surface();
287  assert(tm.recHit()->surface() == &surface);
288  if (!surface.mediumProperties().isValid()) {
289  LogError("TrackingRecoMaterialAnalyser")
290  << "Medium properties for material linked to detector"
291  << " are invalid at: "
292  << current_tsos.globalPosition() << " "
293  << (SiStripDetId)current_det << endl;
294  assert(0);
295  continue;
296  }
297  float p2 = localP.mag2();
298  float xf = std::abs(std::sqrt(p2)/localP.z());
299  float ori_xi = surface.mediumProperties().xi();
300  float ori_radLen = surface.mediumProperties().radLen();
301  float xi = ori_xi*xf;
302  float radLen = ori_radLen*xf;
303 
304  // NOTA BENE: THIS ASSUMES THAT THE TRACKS HAVE BEEN PRODUCED
305  // WITH SPLITTING, AS IS THE CASE FOR THE generalTracks
306  // collection.
307 
308  // Since there are double-sided (glued) modules all over the
309  // tracker, the material budget has been internally
310  // partitioned in two equal components, so that each single
311  // layer will receive half of the correct radLen. For this
312  // reason, only for the double-sided components, we rescale
313  // the obtained radLen by 2.
314 
315  // In particular see code here:
316  // http://cmslxr.fnal.gov/dxr/CMSSW_8_0_5/source/Geometry/TrackerGeometryBuilder/src/TrackerGeomBuilderFromGeometricDet.cc#213
317  // where, in the SiStrip Tracker, if the module has a partner
318  // (i.e. it's a glued detector) the plane is built with a
319  // scaling of 0.5. The actual plane is built few lines below:
320  // http://cmslxr.fnal.gov/dxr/CMSSW_8_0_5/source/Geometry/TrackerGeometryBuilder/src/TrackerGeomBuilderFromGeometricDet.cc#287
321 
322  if (isDoubleSided(current_det)) {
323  LogTrace("TrackingRecoMaterialAnalyser") << "Eta: " << track.eta() << " "
324  << sDETS[current_det.subdetId()] << trk_topology->layer(current_det)
325  << " has ori_radLen: " << ori_radLen << " and ori_xi: " << xi
326  << " and has radLen: " << radLen << " and xi: " << xi << endl;
327  ori_radLen *= 2.;
328  radLen *= 2.;
329  }
330 
331  histosOriEta_[sDETS[current_det.subdetId()]+to_string(trk_topology->layer(current_det))]->Fill(current_tsos.globalPosition().eta(), ori_radLen);
332  histosEta_[sDETS[current_det.subdetId()]+to_string(trk_topology->layer(current_det))]->Fill(current_tsos.globalPosition().eta(), radLen);
333  histo_RZ_Ori_->Fill(current_tsos.globalPosition().z(), current_tsos.globalPosition().perp(), ori_radLen);
334  histo_RZ_->Fill(current_tsos.globalPosition().z(), current_tsos.globalPosition().perp(), radLen);
335  LogInfo("TrackingRecoMaterialAnalyser") << "Eta: " << track.eta() << " "
336  << sDETS[current_det.subdetId()] << trk_topology->layer(current_det)
337  << " has ori_radLen: " << ori_radLen << " and ori_xi: " << xi
338  << " and has radLen: " << radLen << " and xi: " << xi << endl;
339  }
340  }
341  }
342 }
343 
344 //-------------------------------------------------------------------------
345 // define as a plugin
double p() const
momentum vector magnitude
Definition: TrackBase.h:654
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
const edm::EDGetTokenT< reco::VertexCollection > verticesToken_
float radLen() const
T perp() const
Definition: PV3DBase.h:72
void setServices(const edm::EventSetup &) override
set the services needed by the TrackTransformer
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:113
#define nullptr
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
GlobalPoint globalPosition() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::unordered_map< std::string, MonitorElement * > histosEta_
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
const edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:67
const Point & position() const
position
Definition: Vertex.h:109
static const std::vector< std::string > sDETS
unsigned int numberOfLayers(int subdet) const
double outerZ() const
z coordinate of the outermost hit position
Definition: Track.h:154
void Fill(long long x)
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:995
std::unordered_map< std::string, MonitorElement * > histosOriEta_
LocalVector localMomentum() const
const SurfaceType & surface() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:690
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * bookProfile2D(Args &&...args)
Definition: DQMStore.h:114
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
z coordinate
Definition: Vertex.h:115
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
std::vector< Trajectory > transform(const reco::Track &) const override
Convert a reco::Track into Trajectory.
TrackingRecoMaterialAnalyser(const edm::ParameterSet &)
bool isValid() const
Definition: HandleBase.h:74
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
double ndof() const
Definition: Vertex.h:105
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
SubDetector subDetector() const
Definition: SiStripDetId.h:102
Definition: DetId.h:18
uint32_t glued() const
Definition: SiStripDetId.h:158
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:134
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:72
bool isFake() const
Definition: Vertex.h:72
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:485
float xi() const
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:549
unsigned int layer(const DetId &id) const
T eta() const
Definition: PV3DBase.h:76
fixed size matrix
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:62
T get() const
Definition: EventSetup.h:71
bool isDoubleSided(SiStripDetId strip_id)
const Point & position() const
position
Definition: BeamSpot.h:62
void analyze(const edm::Event &, const edm::EventSetup &) override
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:630
bool isValid() const
const MediumProperties & mediumProperties() const
Definition: Surface.h:112
Definition: event.py:1
Definition: Run.h:45