CMS 3D CMS Logo

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