CMS 3D CMS Logo

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