CMS 3D CMS Logo

AlignmentStats.cc
Go to the documentation of this file.
2 
11 
17 
19 
25 
26 using namespace std;
27 
29  : esTokenTTopo_(esConsumes()),
30  esTokenTkGeo_(esConsumes()),
31  src_(iConfig.getParameter<edm::InputTag>("src")),
32  overlapAM_(iConfig.getParameter<edm::InputTag>("OverlapAssoMap")),
33  keepTrackStats_(iConfig.getParameter<bool>("keepTrackStats")),
34  keepHitPopulation_(iConfig.getParameter<bool>("keepHitStats")),
35  statsTreeName_(iConfig.getParameter<string>("TrkStatsFileName")),
36  hitsTreeName_(iConfig.getParameter<string>("HitStatsFileName")),
37  prescale_(iConfig.getParameter<uint32_t>("TrkStatsPrescale")),
38  trackToken_(consumes<reco::TrackCollection>(src_)),
39  mapToken_(consumes<AliClusterValueMap>(overlapAM_)) {
40  //sanity checks
41 
42  //init
43  outtree_ = nullptr;
44 
45 } //end constructor
46 
49  desc.add<edm::InputTag>("src", edm::InputTag("generalTracks"));
50  desc.add<edm::InputTag>("OverlapAssoMap", edm::InputTag("OverlapAssoMap"));
51  desc.add<bool>("keepTrackStats", false);
52  desc.add<bool>("keepHitStats", false);
53  desc.add<std::string>("TrkStatsFileName", "tracks_statistics.root");
54  desc.add<std::string>("HitStatsFileName", "HitMaps.root");
55  desc.add<unsigned int>("TrkStatsPrescale", 1);
56  descriptions.add("AlignmentStats", desc);
57 }
58 
59 void AlignmentStats::beginJob() { // const edm::EventSetup &iSetup
60 
61  //book track stats tree
62  treefile_ = new TFile(statsTreeName_.c_str(), "RECREATE");
63  treefile_->cd();
64  outtree_ = new TTree("AlignmentTrackStats", "Statistics of Tracks used for Alignment");
65  // int nHitsinPXB[MAXTRKS_], nHitsinPXE[MAXTRKS_], nHitsinTEC[MAXTRKS_], nHitsinTIB[MAXTRKS_],nHitsinTOB[MAXTRKS_],nHitsinTID[MAXTRKS_];
66 
67  outtree_->Branch("Ntracks", &ntracks, "Ntracks/i");
68  outtree_->Branch("Run_", &run_, "Run_Nr/I");
69  outtree_->Branch("Event", &event_, "EventNr/I");
70  outtree_->Branch("Eta", &Eta, "Eta[Ntracks]/F");
71  outtree_->Branch("Phi", &Phi, "Phi[Ntracks]/F");
72  outtree_->Branch("P", &P, "P[Ntracks]/F");
73  outtree_->Branch("Pt", &Pt, "Pt[Ntracks]/F");
74  outtree_->Branch("Chi2n", &Chi2n, "Chi2n[Ntracks]/F");
75  outtree_->Branch("Nhits", &Nhits, "Nhits[Ntracks][7]/I");
76  /*
77  outtree_->Branch("NhitsPXB" , ,);
78  outtree_->Branch("NhitsPXE" , ,);
79  outtree_->Branch("NhitsTIB" , ,);
80  outtree_->Branch("NhitsTID" , ,);
81  outtree_->Branch("NhitsTOB" , ,);
82  outtree_->Branch("NhitsTOB" , ,);
83  */
84 
86 } //end beginJob
87 
89  //load list of detunits needed then in endJob
90  if (!trackerGeometry_) {
91  trackerGeometry_ = std::make_unique<TrackerGeometry>(iSetup.getData(esTokenTkGeo_));
92  }
93 
94  if (!trackerTopology_) {
95  trackerTopology_ = std::make_unique<TrackerTopology>(iSetup.getData(esTokenTTopo_));
96  }
97 
98  //take trajectories and tracks to loop on
99  // edm::Handle<TrajTrackAssociationCollection> TrackAssoMap;
101 
102  //take overlap HitAssomap
103  const edm::Handle<AliClusterValueMap> &hMap = iEvent.getHandle(mapToken_);
104  const AliClusterValueMap &OverlapMap = *hMap;
105 
106  // Initialise
107  run_ = 1;
108  event_ = 1;
109  ntracks = 0;
110  run_ = iEvent.id().run();
111  event_ = iEvent.id().event();
112  ntracks = Tracks->size();
113  if (ntracks > 1)
114  edm::LogVerbatim("AlignmenStats") << "~~~~~~~~~~~~\n For this event processing " << ntracks << " tracks";
115 
116  unsigned int trk_cnt = 0;
117 
118  for (int j = 0; j < MAXTRKS_; j++) {
119  Eta[j] = -9999.0;
120  Phi[j] = -8888.0;
121  P[j] = -7777.0;
122  Pt[j] = -6666.0;
123  Chi2n[j] = -2222.0;
124  for (int k = 0; k < 7; k++) {
125  Nhits[j][k] = 0;
126  }
127  }
128 
129  // int npxbhits=0;
130 
131  //loop on tracks
132  for (std::vector<reco::Track>::const_iterator ittrk = Tracks->begin(), edtrk = Tracks->end(); ittrk != edtrk;
133  ++ittrk) {
134  Eta[trk_cnt] = ittrk->eta();
135  Phi[trk_cnt] = ittrk->phi();
136  Chi2n[trk_cnt] = ittrk->normalizedChi2();
137  P[trk_cnt] = ittrk->p();
138  Pt[trk_cnt] = ittrk->pt();
139  Nhits[trk_cnt][0] = ittrk->numberOfValidHits();
140 
141  if (ntracks > 1)
142  edm::LogVerbatim("AlignmenStats") << "Track #" << trk_cnt + 1 << " params: Eta=" << Eta[trk_cnt]
143  << " Phi=" << Phi[trk_cnt] << " P=" << P[trk_cnt]
144  << " Nhits=" << Nhits[trk_cnt][0];
145 
146  int nhit = 0;
147  //loop on tracking rechits
148  //edm::LogVerbatim("AlignmenStats") << " loop on hits of track #" << (itt - tracks->begin());
149  for (auto const &hit : ittrk->recHits()) {
150  if (!hit->isValid())
151  continue;
152  DetId detid = hit->geographicalId();
153  int subDet = detid.subdetId();
154  uint32_t rawId = hit->geographicalId().rawId();
155 
156  // if(subDet==1)npxbhits++;
157 
158  //look if you find this detid in the map
159  DetHitMap::iterator mapiter;
160  mapiter = hitmap_.find(rawId);
161  if (mapiter != hitmap_.end()) { //present, increase its value by one
162  // hitmap_[rawId]=hitmap_[rawId]+1;
163  ++(hitmap_[rawId]);
164  } else { //not present, let's add this key to the map with value=1
165  hitmap_.insert(pair<uint32_t, uint32_t>(rawId, 1));
166  }
167 
168  AlignmentClusterFlag inval;
169 
170  bool hitInPixel = (subDet == PixelSubdetector::PixelBarrel) || (subDet == PixelSubdetector::PixelEndcap);
171  bool hitInStrip = (subDet == SiStripDetId::TIB) || (subDet == SiStripDetId::TID) ||
172  (subDet == SiStripDetId::TOB) || (subDet == SiStripDetId::TEC);
173 
174  if (!(hitInPixel || hitInStrip)) {
175  //skip only this hit, don't stop everything throwing an exception
176  edm::LogError("AlignmentStats") << "Hit not belonging neither to pixels nor to strips! Skipping it. SubDet= "
177  << subDet;
178  continue;
179  }
180 
181  //check also if the hit is an overlap. If yes fill a dedicated hitmap
182  if (hitInStrip) {
183  const std::type_info &type = typeid(*hit);
184 
185  if (type == typeid(SiStripRecHit1D)) {
186  //Notice the difference respect to when one loops on Trajectories: the recHit is a TrackingRecHit and not a TransientTrackingRecHit
187  const SiStripRecHit1D *striphit = dynamic_cast<const SiStripRecHit1D *>(hit);
188  if (striphit != nullptr) {
189  SiStripRecHit1D::ClusterRef stripclust(striphit->cluster());
190  inval = OverlapMap[stripclust];
191  } else {
192  // edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit1D failed! TypeId of the RecHit: "<<className(*hit);
193  throw cms::Exception("NullPointerError")
194  << "ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit1D failed! TypeId of the RecHit: "
195  << className(*hit);
196  }
197  } //end if sistriprechit1D
198  if (type == typeid(SiStripRecHit2D)) {
199  const SiStripRecHit2D *striphit = dynamic_cast<const SiStripRecHit2D *>(hit);
200  if (striphit != nullptr) {
201  SiStripRecHit2D::ClusterRef stripclust(striphit->cluster());
202  inval = OverlapMap[stripclust];
203  //edm::LogVerbatim("AlignmenStats")<<"Taken the Strip Cluster with ProdId "<<stripclust.id() <<"; the Value in the map is "<<inval<<" (DetId is "<<hit->geographicalId().rawId()<<")";
204  } else {
205  throw cms::Exception("NullPointerError")
206  << "ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit2D failed! TypeId of the RecHit: "
207  << className(*hit);
208  // edm::LogError("AlignmentStats")<<"ERROR in <AlignmentStats::analyze>: Dynamic cast of Strip RecHit2D failed! TypeId of the RecHit: "<<className(*hit);
209  }
210  } //end if sistriprechit2D
211 
212  } //end if hit in Strips
213  else {
214  const SiPixelRecHit *pixelhit = dynamic_cast<const SiPixelRecHit *>(hit);
215  if (pixelhit != nullptr) {
216  SiPixelClusterRefNew pixclust(pixelhit->cluster());
217  inval = OverlapMap[pixclust];
218  //edm::LogVerbatim("AlignmenStats")<<"Taken the Pixel Cluster with ProdId "<<pixclust.id() <<"; the Value in the map is "<<inval<<" (DetId is "<<hit->geographicalId().rawId()<<")";
219  } else {
220  edm::LogError("AlignmentStats")
221  << "ERROR in <AlignmentStats::analyze>: Dynamic cast of Pixel RecHit failed! TypeId of the RecHit: "
222  << className(*hit);
223  }
224  } //end else hit is in Pixel
225 
226  bool isOverlapHit(inval.isOverlap());
227 
228  if (isOverlapHit) {
229  edm::LogVerbatim("AlignmenStats") << "This hit is an overlap !";
230  DetHitMap::iterator overlapiter;
231  overlapiter = overlapmap_.find(rawId);
232 
233  if (overlapiter !=
234  overlapmap_.end()) { //the det already collected at least an overlap, increase its value by one
235  overlapmap_[rawId] = overlapmap_[rawId] + 1;
236  } else { //first overlap on det unit, let's add it to the map
237  overlapmap_.insert(pair<uint32_t, uint32_t>(rawId, 1));
238  }
239  } //end if the hit is an overlap
240 
241  int subdethit = static_cast<int>(hit->geographicalId().subdetId());
242  if (ntracks > 1)
243  edm::LogVerbatim("AlignmenStats") << "Hit in SubDet=" << subdethit;
244  Nhits[trk_cnt][subdethit] = Nhits[trk_cnt][subdethit] + 1;
245  nhit++;
246  } //end loop on trackingrechits
247  trk_cnt++;
248 
249  } //end loop on tracks
250 
251  //edm::LogVerbatim("AlignmenStats")<<"Total number of pixel hits is " << npxbhits;
252 
253  tmpPresc_--;
254  if (tmpPresc_ < 1) {
255  outtree_->Fill();
257  }
258  if (trk_cnt != ntracks)
259  edm::LogError("AlignmentStats") << "\nERROR! trk_cnt=" << trk_cnt << " ntracks=" << ntracks;
260 
261  return;
262 }
263 
265  treefile_->cd();
266  edm::LogInfo("AlignmentStats") << "Writing out the TrackStatistics in " << gDirectory->GetPath();
267  outtree_->Write();
268  delete outtree_;
269 
270  //create tree with hit maps (hitstree)
271  //book track stats tree
272  TFile *hitsfile = new TFile(hitsTreeName_.c_str(), "RECREATE");
273  hitsfile->cd();
274  TTree *hitstree = new TTree("AlignmentHitMap", "Maps of Hits used for Alignment");
275 
276  unsigned int id = 0, nhits = 0, noverlaps = 0;
277  float posX(-99999.0), posY(-77777.0), posZ(-88888.0);
278  float posEta(-6666.0), posPhi(-5555.0), posR(-4444.0);
279  int subdet = 0;
280  unsigned int layer = 0;
281  bool is2D = false, isStereo = false;
282  hitstree->Branch("DetId", &id, "DetId/i");
283  hitstree->Branch("Nhits", &nhits, "Nhits/i");
284  hitstree->Branch("Noverlaps", &noverlaps, "Noverlaps/i");
285  hitstree->Branch("SubDet", &subdet, "SubDet/I");
286  hitstree->Branch("Layer", &layer, "Layer/i");
287  hitstree->Branch("is2D", &is2D, "is2D/B");
288  hitstree->Branch("isStereo", &isStereo, "isStereo/B");
289  hitstree->Branch("posX", &posX, "posX/F");
290  hitstree->Branch("posY", &posY, "posY/F");
291  hitstree->Branch("posZ", &posZ, "posZ/F");
292  hitstree->Branch("posR", &posR, "posR/F");
293  hitstree->Branch("posEta", &posEta, "posEta/F");
294  hitstree->Branch("posPhi", &posPhi, "posPhi/F");
295 
296  /*
297  TTree *overlapstree=new TTree("OverlapHitMap","Maps of Overlaps used for Alignment");
298  hitstree->Branch("DetId", &id , "DetId/i");
299  hitstree->Branch("NOverlaps", &nhits , "Nhits/i");
300  hitstree->Branch("SubDet", &subdet, "SubDet/I");
301  hitstree->Branch("Layer", &layer, "Layer/i");
302  hitstree->Branch("is2D" , &is2D, "is2D/B");
303  hitstree->Branch("isStereo",&isStereo,"isStereo/B");
304  hitstree->Branch("posX", &posX, "posX/F");
305  hitstree->Branch("posY", &posY, "posY/F");
306  hitstree->Branch("posZ", &posZ, "posZ/F");
307  hitstree->Branch("posR", &posR, "posR/F");
308  hitstree->Branch("posEta", &posEta, "posEta/F");
309  hitstree->Branch("posPhi", &posPhi, "posPhi/F");
310  */
311 
312  std::unique_ptr<AlignableTracker> theAliTracker =
313  std::make_unique<AlignableTracker>(trackerGeometry_.get(), trackerTopology_.get());
314  const auto &Detunitslist = theAliTracker->deepComponents();
315  int ndetunits = Detunitslist.size();
316  edm::LogInfo("AlignmentStats") << "Number of DetUnits in the AlignableTracker: " << ndetunits;
317 
318  for (int det_cnt = 0; det_cnt < ndetunits; ++det_cnt) {
319  //re-initialize for safety
320  id = 0;
321  nhits = 0;
322  noverlaps = 0;
323  posX = -99999.0;
324  posY = -77777.0;
325  posZ = -88888.0;
326  posEta = -6666.0;
327  posPhi = -5555.0;
328  posR = -4444.0;
329  subdet = 0;
330  layer = 0;
331  is2D = false;
332  isStereo = false;
333 
334  //if detunit in vector is found also in the map, look for how many hits were collected
335  //and save in the tree this number
336  id = static_cast<uint32_t>(Detunitslist[det_cnt]->id());
337  if (hitmap_.find(id) != hitmap_.end()) {
338  nhits = hitmap_[id];
339  }
340  //if not, save nhits=0
341  else {
342  nhits = 0;
343  hitmap_.insert(pair<uint32_t, uint32_t>(id, 0));
344  }
345 
346  if (overlapmap_.find(id) != overlapmap_.end()) {
347  noverlaps = overlapmap_[id];
348  }
349  //if not, save nhits=0
350  else {
351  noverlaps = 0;
352  overlapmap_.insert(pair<uint32_t, uint32_t>(id, 0));
353  }
354 
355  //take other geometrical infos from the det
356  posX = Detunitslist[det_cnt]->globalPosition().x();
357  posY = Detunitslist[det_cnt]->globalPosition().y();
358  posZ = Detunitslist[det_cnt]->globalPosition().z();
359 
360  align::GlobalVector vec(posX, posY, posZ);
361  posR = vec.perp();
362  posPhi = vec.phi();
363  posEta = vec.eta();
364  // posPhi = atan2(posY,posX);
365 
366  DetId detid(id);
367  subdet = detid.subdetId();
368 
369  //get layers, petals, etc...
370  if (subdet == PixelSubdetector::PixelBarrel) { //PXB
371 
372  layer = trackerTopology_->pxbLayer(id);
373  is2D = true;
374  isStereo = false;
375  } else if (subdet == PixelSubdetector::PixelEndcap) {
376  layer = trackerTopology_->pxfDisk(id);
377  is2D = true;
378  isStereo = false;
379  } else if (subdet == SiStripDetId::TIB) {
380  layer = trackerTopology_->tibLayer(id);
381  is2D = trackerTopology_->tibIsDoubleSide(id);
382  isStereo = trackerTopology_->tibIsStereo(id);
383  } else if (subdet == SiStripDetId::TID) {
384  layer = trackerTopology_->tidWheel(id);
385  is2D = trackerTopology_->tidIsDoubleSide(id);
386  isStereo = trackerTopology_->tidIsStereo(id);
387  } else if (subdet == SiStripDetId::TOB) {
388  layer = trackerTopology_->tobLayer(id);
389  is2D = trackerTopology_->tobIsDoubleSide(id);
390  isStereo = trackerTopology_->tobIsStereo(id);
391  } else if (subdet == SiStripDetId::TEC) {
392  layer = trackerTopology_->tecWheel(id);
393  is2D = trackerTopology_->tecIsDoubleSide(id);
394  isStereo = trackerTopology_->tecIsStereo(id);
395  } else {
396  edm::LogError("AlignmentStats") << "Detector not belonging neither to pixels nor to strips! Skipping it. SubDet= "
397  << subdet;
398  }
399 
400  //write in the hitstree
401  hitstree->Fill();
402  } //end loop over detunits
403 
404  //save hitstree
405  hitstree->Write();
406  delete hitstree;
407  //delete Detunitslist;
408  hitmap_.clear();
409  overlapmap_.clear();
410  delete hitsfile;
411 }
412 // ========= MODULE DEF ==============
float Eta[MAXTRKS_]
Log< level::Info, true > LogVerbatim
DetHitMap hitmap_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const std::string hitsTreeName_
T perp() const
Definition: PV3DBase.h:69
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > esTokenTkGeo_
static constexpr auto TID
Definition: SiStripDetId.h:38
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
DetHitMap overlapmap_
T eta() const
Definition: PV3DBase.h:73
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
float Chi2n[MAXTRKS_]
const uint32_t prescale_
Log< level::Error, false > LogError
ClusterRef cluster() const
std::unique_ptr< TrackerGeometry > trackerGeometry_
const edm::EDGetTokenT< reco::TrackCollection > trackToken_
constexpr std::array< uint8_t, layerIndexSize > layer
int iEvent
Definition: GenABIO.cc:224
uint32_t tmpPresc_
bool is2D(HitType hitType)
float Phi[MAXTRKS_]
#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
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static constexpr auto TOB
Definition: SiStripDetId.h:39
Log< level::Info, false > LogInfo
Definition: DetId.h:17
const std::string statsTreeName_
const edm::EDGetTokenT< AliClusterValueMap > mapToken_
AlignmentStats(const edm::ParameterSet &iConfig)
static const int MAXTRKS_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::pair< OmniClusterRef, TrackingParticleRef > P
static constexpr auto TIB
Definition: SiStripDetId.h:37
unsigned int ntracks
fixed size matrix
int Nhits[MAXTRKS_][7]
HLT enums.
float Pt[MAXTRKS_]
void beginJob() override
ClusterRef cluster() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::unique_ptr< TrackerTopology > trackerTopology_
static constexpr auto TEC
Definition: SiStripDetId.h:40
std::string className(const T &t)
Definition: ClassName.h:31
void endJob() override
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > esTokenTTopo_
Our base class.
Definition: SiPixelRecHit.h:23
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override