CMS 3D CMS Logo

ApeEstimator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ApeEstimator
4 // Class: ApeEstimator
5 //
13 //
14 // Original Author: Johannes Hauk
15 // Created: Tue Jan 6 15:02:09 CET 2009
16 // Modified by: Christian Schomakers (RWTH Aachen)
17 // $Id: ApeEstimator.cc,v 1.27 2012/06/26 09:42:33 hauk Exp $
18 //
19 //
20 
21 // system include files
22 #include <memory>
23 #include <sstream>
24 #include <fstream>
25 
26 // user include files
37 
40 
60 
68 
72 //added by Ajay 6Nov 2014
73 //.......................
76 
78 
79 //...............
80 //
81 
84 
86 
88 
93 
94 #include "TH1.h"
95 #include "TH2.h"
96 #include "TProfile.h"
97 #include "TFile.h"
98 #include "TTree.h"
99 #include "TF1.h"
100 #include "TString.h"
101 #include "TMath.h"
102 
105 
108 //ADDED BY LOIC QUERTENMONT
120 
121 //
122 // class decleration
123 //
124 
126 public:
127  explicit ApeEstimator(const edm::ParameterSet&);
128  ~ApeEstimator() override;
129 
130 private:
132  PositionAndError2() : posX(-999.F), posY(-999.F), errX2(-999.F), errY2(-999.F){};
133  PositionAndError2(float x, float y, float eX, float eY) : posX(x), posY(y), errX2(eX), errY2(eY){};
134  float posX;
135  float posY;
136  float errX2;
137  float errY2;
138  };
139  typedef std::pair<TrackStruct::HitState, PositionAndError2> StatePositionAndError2;
140 
141  void beginJob() override;
142  void analyze(const edm::Event&, const edm::EventSetup&) override;
143  void endJob() override;
144 
145  bool isHit2D(const TrackingRecHit&) const;
146 
147  void sectorBuilder();
148  bool checkIntervalsForSectors(const unsigned int sectorCounter, const std::vector<double>&) const;
149  bool checkModuleIds(const unsigned int, const std::vector<unsigned int>&) const;
150  bool checkModuleBools(const bool, const std::vector<unsigned int>&) const;
151  bool checkModuleDirections(const int, const std::vector<int>&) const;
152  bool checkModulePositions(const float, const std::vector<double>&) const;
153  void statistics(const TrackerSectorStruct&, const int) const;
154 
155  void residualErrorBinning();
156 
159  void bookTrackHists();
160 
163 
164  StatePositionAndError2 positionAndError2(const LocalPoint&, const LocalError&, const TransientTrackingRecHit&);
167 
168  void hitSelection();
169  void setHitSelectionMap(const std::string&);
170  void setHitSelectionMapUInt(const std::string&);
172  bool inDoubleInterval(const std::vector<double>&, const float) const;
173  bool inUintInterval(const std::vector<unsigned int>&, const unsigned int, const unsigned int = 999) const;
174 
179 
180  void calculateAPE();
181 
182  // ----------member data ---------------------------
184  std::map<unsigned int, TrackerSectorStruct> m_tkSector_;
186 
189 
190  std::map<unsigned int, std::pair<double, double> > m_resErrBins_;
191  std::map<unsigned int, ReducedTrackerTreeVariables> m_tkTreeVar_;
192 
193  std::map<std::string, std::vector<double> > m_hitSelection_;
194  std::map<std::string, std::vector<unsigned int> > m_hitSelectionUInt_;
195 
196  bool trackCut_;
197 
198  const unsigned int maxTracksPerEvent_;
199  const unsigned int minGoodHitsPerTrack_;
200 
201  const bool analyzerMode_;
202 
203  const bool calculateApe_;
204 
206 };
207 
208 //
209 // constants, enums and typedefs
210 //
211 
212 //
213 // static data member definitions
214 //
215 
216 //
217 // constructors and destructor
218 //
220  : parameterSet_(iConfig),
221  tjTagToken_(
222  consumes<TrajTrackAssociationCollection>(parameterSet_.getParameter<edm::InputTag>("tjTkAssociationMapTag"))),
223  offlinebeamSpot_(consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"))),
224  trackCut_(false),
225  maxTracksPerEvent_(parameterSet_.getParameter<unsigned int>("maxTracksPerEvent")),
226  minGoodHitsPerTrack_(parameterSet_.getParameter<unsigned int>("minGoodHitsPerTrack")),
227  analyzerMode_(parameterSet_.getParameter<bool>("analyzerMode")),
228  calculateApe_(parameterSet_.getParameter<bool>("calculateApe")) {
230 }
231 
233 
234 //
235 // member functions
236 //
237 
238 // -----------------------------------------------------------------------------------------------------------
239 
241  TFile* tkTreeFile(TFile::Open((parameterSet_.getParameter<std::string>("TrackerTreeFile")).c_str()));
242  if (tkTreeFile) {
243  edm::LogInfo("SectorBuilder") << "TrackerTreeFile OK";
244  } else {
245  edm::LogError("SectorBuilder") << "TrackerTreeFile not found";
246  return;
247  }
248  TTree* tkTree(nullptr);
249  tkTreeFile->GetObject("TrackerTreeGenerator/TrackerTree/TrackerTree", tkTree);
250  if (tkTree) {
251  edm::LogInfo("SectorBuilder") << "TrackerTree OK";
252  } else {
253  edm::LogError("SectorBuilder") << "TrackerTree not found in file";
254  return;
255  }
256  unsigned int rawId(999), subdetId(999), layer(999), side(999), half(999), rod(999), ring(999), petal(999), blade(999),
257  panel(999), outerInner(999), module(999), nStrips(999);
258  bool isDoubleSide(false), isRPhi(false), isStereo(false);
259  int uDirection(999), vDirection(999), wDirection(999);
260  float posR(999.F), posPhi(999.F), posEta(999.F), posX(999.F), posY(999.F), posZ(999.F);
261 
262  tkTree->SetBranchAddress("RawId", &rawId);
263  tkTree->SetBranchAddress("SubdetId", &subdetId);
264  tkTree->SetBranchAddress("Layer", &layer);
265  tkTree->SetBranchAddress("Side", &side);
266  tkTree->SetBranchAddress("Half", &half);
267  tkTree->SetBranchAddress("Rod", &rod);
268  tkTree->SetBranchAddress("Ring", &ring);
269  tkTree->SetBranchAddress("Petal", &petal);
270  tkTree->SetBranchAddress("Blade", &blade);
271  tkTree->SetBranchAddress("Panel", &panel);
272  tkTree->SetBranchAddress("OuterInner", &outerInner);
273  tkTree->SetBranchAddress("Module", &module);
274  tkTree->SetBranchAddress("NStrips", &nStrips);
275  tkTree->SetBranchAddress("IsDoubleSide", &isDoubleSide);
276  tkTree->SetBranchAddress("IsRPhi", &isRPhi);
277  tkTree->SetBranchAddress("IsStereo", &isStereo);
278  tkTree->SetBranchAddress("UDirection", &uDirection);
279  tkTree->SetBranchAddress("VDirection", &vDirection);
280  tkTree->SetBranchAddress("WDirection", &wDirection);
281  tkTree->SetBranchAddress("PosR", &posR);
282  tkTree->SetBranchAddress("PosPhi", &posPhi);
283  tkTree->SetBranchAddress("PosEta", &posEta);
284  tkTree->SetBranchAddress("PosX", &posX);
285  tkTree->SetBranchAddress("PosY", &posY);
286  tkTree->SetBranchAddress("PosZ", &posZ);
287 
288  int nModules(tkTree->GetEntries());
289  TrackerSectorStruct allSectors;
290 
291  //Loop over all Sectors
292  unsigned int sectorCounter(0);
293  std::vector<edm::ParameterSet> v_sectorDef(parameterSet_.getParameter<std::vector<edm::ParameterSet> >("Sectors"));
294  edm::LogInfo("SectorBuilder") << "There are " << v_sectorDef.size() << " Sectors definded";
295  for (auto const& parSet : v_sectorDef) {
296  ++sectorCounter;
297  const std::string& sectorName(parSet.getParameter<std::string>("name"));
298  std::vector<unsigned int> v_rawId(parSet.getParameter<std::vector<unsigned int> >("rawId")),
299  v_subdetId(parSet.getParameter<std::vector<unsigned int> >("subdetId")),
300  v_layer(parSet.getParameter<std::vector<unsigned int> >("layer")),
301  v_side(parSet.getParameter<std::vector<unsigned int> >("side")),
302  v_half(parSet.getParameter<std::vector<unsigned int> >("half")),
303  v_rod(parSet.getParameter<std::vector<unsigned int> >("rod")),
304  v_ring(parSet.getParameter<std::vector<unsigned int> >("ring")),
305  v_petal(parSet.getParameter<std::vector<unsigned int> >("petal")),
306  v_blade(parSet.getParameter<std::vector<unsigned int> >("blade")),
307  v_panel(parSet.getParameter<std::vector<unsigned int> >("panel")),
308  v_outerInner(parSet.getParameter<std::vector<unsigned int> >("outerInner")),
309  v_module(parSet.getParameter<std::vector<unsigned int> >("module")),
310  v_nStrips(parSet.getParameter<std::vector<unsigned int> >("nStrips")),
311  v_isDoubleSide(parSet.getParameter<std::vector<unsigned int> >("isDoubleSide")),
312  v_isRPhi(parSet.getParameter<std::vector<unsigned int> >("isRPhi")),
313  v_isStereo(parSet.getParameter<std::vector<unsigned int> >("isStereo"));
314  std::vector<int> v_uDirection(parSet.getParameter<std::vector<int> >("uDirection")),
315  v_vDirection(parSet.getParameter<std::vector<int> >("vDirection")),
316  v_wDirection(parSet.getParameter<std::vector<int> >("wDirection"));
317  std::vector<double> v_posR(parSet.getParameter<std::vector<double> >("posR")),
318  v_posPhi(parSet.getParameter<std::vector<double> >("posPhi")),
319  v_posEta(parSet.getParameter<std::vector<double> >("posEta")),
320  v_posX(parSet.getParameter<std::vector<double> >("posX")),
321  v_posY(parSet.getParameter<std::vector<double> >("posY")),
322  v_posZ(parSet.getParameter<std::vector<double> >("posZ"));
323 
324  if (!this->checkIntervalsForSectors(sectorCounter, v_posR) ||
325  !this->checkIntervalsForSectors(sectorCounter, v_posPhi) ||
326  !this->checkIntervalsForSectors(sectorCounter, v_posEta) ||
327  !this->checkIntervalsForSectors(sectorCounter, v_posX) ||
328  !this->checkIntervalsForSectors(sectorCounter, v_posY) ||
329  !this->checkIntervalsForSectors(sectorCounter, v_posZ))
330  continue;
331 
332  TrackerSectorStruct tkSector;
333  tkSector.name = sectorName;
334 
335  ReducedTrackerTreeVariables tkTreeVar;
336 
337  //Loop over all Modules
338  for (int module = 0; module < nModules; ++module) {
339  tkTree->GetEntry(module);
340 
341  if (sectorCounter == 1) {
342  tkTreeVar.subdetId = subdetId;
343  tkTreeVar.nStrips = nStrips;
344  tkTreeVar.uDirection = uDirection;
345  tkTreeVar.vDirection = vDirection;
346  tkTreeVar.wDirection = wDirection;
347  m_tkTreeVar_[rawId] = tkTreeVar;
348  }
349 
350  if (!this->checkModuleIds(rawId, v_rawId))
351  continue;
352  if (!this->checkModuleIds(subdetId, v_subdetId))
353  continue;
354  if (!this->checkModuleIds(layer, v_layer))
355  continue;
356  if (!this->checkModuleIds(side, v_side))
357  continue;
358  if (!this->checkModuleIds(half, v_half))
359  continue;
360  if (!this->checkModuleIds(rod, v_rod))
361  continue;
362  if (!this->checkModuleIds(ring, v_ring))
363  continue;
364  if (!this->checkModuleIds(petal, v_petal))
365  continue;
366  if (!this->checkModuleIds(blade, v_blade))
367  continue;
368  if (!this->checkModuleIds(panel, v_panel))
369  continue;
370  if (!this->checkModuleIds(outerInner, v_outerInner))
371  continue;
372  if (!this->checkModuleIds(module, v_module))
373  continue;
374  if (!this->checkModuleIds(nStrips, v_nStrips))
375  continue;
376  if (!this->checkModuleBools(isDoubleSide, v_isDoubleSide))
377  continue;
378  if (!this->checkModuleBools(isRPhi, v_isRPhi))
379  continue;
380  if (!this->checkModuleBools(isStereo, v_isStereo))
381  continue;
382  if (!this->checkModuleDirections(uDirection, v_uDirection))
383  continue;
384  if (!this->checkModuleDirections(vDirection, v_vDirection))
385  continue;
386  if (!this->checkModuleDirections(wDirection, v_wDirection))
387  continue;
388  if (!this->checkModulePositions(posR, v_posR))
389  continue;
390  if (!this->checkModulePositions(posPhi, v_posPhi))
391  continue;
392  if (!this->checkModulePositions(posEta, v_posEta))
393  continue;
394  if (!this->checkModulePositions(posX, v_posX))
395  continue;
396  if (!this->checkModulePositions(posY, v_posY))
397  continue;
398  if (!this->checkModulePositions(posZ, v_posZ))
399  continue;
400 
401  tkSector.v_rawId.push_back(rawId);
402  bool moduleSelected(false);
403  for (auto const& i_rawId : allSectors.v_rawId) {
404  if (rawId == i_rawId)
405  moduleSelected = true;
406  }
407  if (!moduleSelected)
408  allSectors.v_rawId.push_back(rawId);
409  }
410 
411  bool isPixel(false);
412  bool isStrip(false);
413  for (auto const& i_rawId : tkSector.v_rawId) {
414  switch (m_tkTreeVar_[i_rawId].subdetId) {
417  isPixel = true;
418  break;
423  isStrip = true;
424  break;
425  }
426  }
427 
428  if (isPixel && isStrip) {
429  edm::LogError("SectorBuilder")
430  << "Incorrect Sector Definition: there are pixel and strip modules within one sector"
431  << "\n... sector selection is not applied, sector " << sectorCounter << " is not built";
432  continue;
433  }
434  tkSector.isPixel = isPixel;
435 
436  m_tkSector_[sectorCounter] = tkSector;
437  edm::LogInfo("SectorBuilder") << "There are " << tkSector.v_rawId.size() << " Modules in Sector " << sectorCounter;
438  }
439  this->statistics(allSectors, nModules);
440  return;
441 }
442 
443 // -----------------------------------------------------------------------------------------------------------
444 
445 bool ApeEstimator::checkIntervalsForSectors(const unsigned int sectorCounter, const std::vector<double>& v_id) const {
446  if (v_id.empty())
447  return true;
448  if (v_id.size() % 2 == 1) {
449  edm::LogError("SectorBuilder")
450  << "Incorrect Sector Definition: Position Vectors need even number of arguments (Intervals)"
451  << "\n... sector selection is not applied, sector " << sectorCounter << " is not built";
452  return false;
453  }
454  int entry(0);
455  double intervalBegin(999.);
456  for (auto const& i_id : v_id) {
457  ++entry;
458  if (entry % 2 == 1)
459  intervalBegin = i_id;
460  if (entry % 2 == 0 && intervalBegin > i_id) {
461  edm::LogError("SectorBuilder") << "Incorrect Sector Definition (Position Vector Intervals): \t" << intervalBegin
462  << " is bigger than " << i_id << " but is expected to be smaller"
463  << "\n... sector selection is not applied, sector " << sectorCounter
464  << " is not built";
465  return false;
466  }
467  }
468  return true;
469 }
470 
471 bool ApeEstimator::checkModuleIds(const unsigned int id, const std::vector<unsigned int>& v_id) const {
472  if (v_id.empty())
473  return true;
474  for (auto const& i_id : v_id) {
475  if (id == i_id)
476  return true;
477  }
478  return false;
479 }
480 
481 bool ApeEstimator::checkModuleBools(const bool id, const std::vector<unsigned int>& v_id) const {
482  if (v_id.empty())
483  return true;
484  for (auto const& i_id : v_id) {
485  if (1 == i_id && id)
486  return true;
487  if (2 == i_id && !id)
488  return true;
489  }
490  return false;
491 }
492 
493 bool ApeEstimator::checkModuleDirections(const int id, const std::vector<int>& v_id) const {
494  if (v_id.empty())
495  return true;
496  for (auto const& i_id : v_id) {
497  if (id == i_id)
498  return true;
499  }
500  return false;
501 }
502 
503 bool ApeEstimator::checkModulePositions(const float id, const std::vector<double>& v_id) const {
504  if (v_id.empty())
505  return true;
506  int entry(0);
507  double intervalBegin(999.);
508  for (auto const& i_id : v_id) {
509  ++entry;
510  if (entry % 2 == 1)
511  intervalBegin = i_id;
512  if (entry % 2 == 0 && id >= intervalBegin && id < i_id)
513  return true;
514  }
515  return false;
516 }
517 
518 void ApeEstimator::statistics(const TrackerSectorStruct& allSectors, const int nModules) const {
519  bool commonModules(false);
520  for (std::map<unsigned int, TrackerSectorStruct>::const_iterator i_sector = m_tkSector_.begin();
521  i_sector != m_tkSector_.end();
522  ++i_sector) {
523  std::map<unsigned int, TrackerSectorStruct>::const_iterator i_sector2(i_sector);
524  for (++i_sector2; i_sector2 != m_tkSector_.end(); ++i_sector2) {
525  unsigned int nCommonModules(0);
526  for (auto const& i_module : (*i_sector).second.v_rawId) {
527  for (auto const& i_module2 : (*i_sector2).second.v_rawId) {
528  if (i_module2 == i_module)
529  ++nCommonModules;
530  }
531  }
532  if (nCommonModules == 0)
533  ; //edm::LogInfo("SectorBuilder")<<"Sector "<<(*i_sector).first<<" and Sector "<<(*i_sector2).first<< " have ZERO Modules in common";
534  else {
535  edm::LogError("SectorBuilder") << "Sector " << (*i_sector).first << " and Sector " << (*i_sector2).first
536  << " have " << nCommonModules << " Modules in common";
537  commonModules = true;
538  }
539  }
540  }
541  if (static_cast<int>(allSectors.v_rawId.size()) == nModules)
542  edm::LogInfo("SectorBuilder") << "ALL Tracker Modules are contained in the Sectors";
543  else
544  edm::LogWarning("SectorBuilder") << "There are " << allSectors.v_rawId.size() << " Modules in all Sectors"
545  << " out of " << nModules << " Tracker Modules";
546  if (!commonModules)
547  edm::LogInfo("SectorBuilder") << "There are ZERO modules associated to different sectors, no ambiguities exist";
548  else
549  edm::LogError("SectorBuilder")
550  << "There are modules associated to different sectors, APE value cannot be assigned reasonably";
551 }
552 
553 // -----------------------------------------------------------------------------------------------------------
554 
556  std::vector<double> v_residualErrorBinning(parameterSet_.getParameter<std::vector<double> >("residualErrorBinning"));
557  if (v_residualErrorBinning.size() == 1) {
558  edm::LogError("ResidualErrorBinning") << "Incorrect selection of Residual Error Bins (used for APE calculation): \t"
559  << "Only one argument passed, so no interval is specified"
560  << "\n... delete whole bin selection"; //m_resErrBins_ remains empty
561  return;
562  }
563  double xMin(0.), xMax(0.);
564  unsigned int binCounter(-1);
565  for (auto const& i_binning : v_residualErrorBinning) {
566  ++binCounter;
567  if (binCounter == 0) {
568  xMin = i_binning;
569  continue;
570  }
571  xMax = i_binning;
572  if (xMax <= xMin) {
573  edm::LogError("ResidualErrorBinning")
574  << "Incorrect selection of Residual Error Bins (used for APE calculation): \t" << xMin << " is bigger than "
575  << xMax << " but is expected to be smaller"
576  << "\n... delete whole bin selection";
577  m_resErrBins_.clear();
578  return;
579  }
580  m_resErrBins_[binCounter].first = xMin;
581  m_resErrBins_[binCounter].second = xMax;
582  xMin = xMax;
583  }
584  edm::LogInfo("ResidualErrorBinning")
585  << m_resErrBins_.size() << " Intervals of residual errors used for separate APE calculation sucessfully set";
586 }
587 
588 // -----------------------------------------------------------------------------------------------------------
589 
591  std::vector<unsigned int> v_errHists(parameterSet_.getParameter<std::vector<unsigned int> >("vErrHists"));
592  for (std::vector<unsigned int>::iterator i_errHists = v_errHists.begin(); i_errHists != v_errHists.end();
593  ++i_errHists) {
594  for (std::vector<unsigned int>::iterator i_errHists2 = i_errHists; i_errHists2 != v_errHists.end();) {
595  ++i_errHists2;
596  if (*i_errHists == *i_errHists2) {
597  edm::LogError("BookSectorHists") << "Value of vErrHists in config exists twice: " << *i_errHists
598  << "\n... delete one of both";
599  v_errHists.erase(i_errHists2);
600  }
601  }
602  }
603 
604  for (auto& i_sector : m_tkSector_) {
605  bool zoomHists(parameterSet_.getParameter<bool>("zoomHists"));
606 
607  double widthMax = zoomHists ? 20. : 200.;
608  double chargePixelMax = zoomHists ? 200000. : 2000000.;
609  double chargeStripMax = zoomHists ? 1000. : 10000.;
610  double sOverNMax = zoomHists ? 200. : 2000.;
611  double logClusterProbMin = zoomHists ? -5. : -15.;
612 
613  double resXAbsMax = zoomHists ? 0.5 : 5.;
614  double norResXAbsMax = zoomHists ? 10. : 50.;
615  double probXMin = zoomHists ? -0.01 : -0.1;
616  double probXMax = zoomHists ? 0.11 : 1.1;
617  double sigmaXMin = zoomHists ? 0. : -0.05;
618  double sigmaXMax = zoomHists ? 0.02 : 1.;
619  double sigmaX2Max = sigmaXMax * sigmaXMax;
620  double sigmaXHitMax = zoomHists ? 0.02 : 1.;
621  double phiSensXMax = zoomHists ? 31. : 93.;
622 
623  double norChi2Max = zoomHists ? 5. : 1000.;
624  double d0Max = zoomHists ? 0.02 : 40.; // cosmics: 100.|100.
625  double dzMax = zoomHists ? 15. : 100.; // cosmics: 200.|600.
626  double pMax = zoomHists ? 200. : 2000.;
627  double invPMax = zoomHists ? 0.05 : 10.; //begins at 20GeV, 0.1GeV
628 
630  if (!fileService) {
631  throw edm::Exception(edm::errors::Configuration, "TFileService is not registered in cfg file");
632  }
633 
634  std::stringstream sector;
635  sector << "Sector_" << i_sector.first;
636  TFileDirectory secDir = fileService->mkdir(sector.str());
637 
638  // Dummy histo containing the sector name as title
639  i_sector.second.Name = secDir.make<TH1F>("z_name", i_sector.second.name.c_str(), 1, 0, 1);
640 
641  // Do not book histos for empty sectors
642  if (i_sector.second.v_rawId.empty()) {
643  continue;
644  }
645  // Set parameters for correlationHists
646  i_sector.second.setCorrHistParams(&secDir, norResXAbsMax, sigmaXHitMax, sigmaXMax);
647 
648  // Book pixel or strip specific hists
649  const bool pixelSector(i_sector.second.isPixel);
650 
651  // Cluster Parameters
652  i_sector.second.m_correlationHistsX["WidthX"] = i_sector.second.bookCorrHistsX("WidthX",
653  "cluster width",
654  "w_{cl,x}",
655  "[# channels]",
656  static_cast<int>(widthMax),
657  static_cast<int>(widthMax),
658  0.,
659  widthMax,
660  "nph");
661  i_sector.second.m_correlationHistsX["BaryStripX"] = i_sector.second.bookCorrHistsX(
662  "BaryStripX", "barycenter of cluster charge", "b_{cl,x}", "[# channels]", 800, 100, -10., 790., "nph");
663 
664  if (pixelSector) {
665  i_sector.second.m_correlationHistsY["WidthY"] = i_sector.second.bookCorrHistsY("WidthY",
666  "cluster width",
667  "w_{cl,y}",
668  "[# channels]",
669  static_cast<int>(widthMax),
670  static_cast<int>(widthMax),
671  0.,
672  widthMax,
673  "nph");
674  i_sector.second.m_correlationHistsY["BaryStripY"] = i_sector.second.bookCorrHistsY(
675  "BaryStripY", "barycenter of cluster charge", "b_{cl,y}", "[# channels]", 800, 100, -10., 790., "nph");
676 
677  i_sector.second.m_correlationHistsX["ChargePixel"] = i_sector.second.bookCorrHistsX(
678  "ChargePixel", "cluster charge", "c_{cl}", "[e]", 100, 50, 0., chargePixelMax, "nph");
679  i_sector.second.m_correlationHistsX["ClusterProbXY"] = i_sector.second.bookCorrHistsX(
680  "ClusterProbXY", "cluster probability xy", "prob_{cl,xy}", "", 100, 50, 0., 1., "nph");
681  i_sector.second.m_correlationHistsX["ClusterProbQ"] = i_sector.second.bookCorrHistsX(
682  "ClusterProbQ", "cluster probability q", "prob_{cl,q}", "", 100, 50, 0., 1., "nph");
683  i_sector.second.m_correlationHistsX["ClusterProbXYQ"] = i_sector.second.bookCorrHistsX(
684  "ClusterProbXYQ", "cluster probability xyq", "prob_{cl,xyq}", "", 100, 50, 0., 1., "nph");
685  i_sector.second.m_correlationHistsX["LogClusterProb"] = i_sector.second.bookCorrHistsX(
686  "LogClusterProb", "cluster probability xy", "log(prob_{cl,xy})", "", 60, 30, logClusterProbMin, 0., "nph");
687  i_sector.second.m_correlationHistsX["IsOnEdge"] =
688  i_sector.second.bookCorrHistsX("IsOnEdge", "IsOnEdge", "isOnEdge", "", 2, 2, 0, 2, "nph");
689  i_sector.second.m_correlationHistsX["HasBadPixels"] =
690  i_sector.second.bookCorrHistsX("HasBadPixels", "HasBadPixels", "hasBadPixels", "", 2, 2, 0, 2, "nph");
691  i_sector.second.m_correlationHistsX["SpansTwoRoc"] =
692  i_sector.second.bookCorrHistsX("SpansTwoRoc", "SpansTwoRoc", "spansTwoRoc", "", 2, 2, 0, 2, "nph");
693  i_sector.second.m_correlationHistsX["QBin"] =
694  i_sector.second.bookCorrHistsX("QBin", "q bin", "q bin", "", 8, 8, 0, 8, "nph");
695 
696  i_sector.second.m_correlationHistsY["ChargePixel"] = i_sector.second.bookCorrHistsY(
697  "ChargePixel", "cluster charge", "c_{cl}", "[e]", 100, 50, 0., chargePixelMax, "nph");
698  i_sector.second.m_correlationHistsY["ClusterProbXY"] = i_sector.second.bookCorrHistsY(
699  "ClusterProbXY", "cluster probability xy", "prob_{cl,xy}", "", 100, 50, 0., 1., "nph");
700  i_sector.second.m_correlationHistsY["ClusterProbQ"] = i_sector.second.bookCorrHistsY(
701  "ClusterProbQ", "cluster probability q", "prob_{cl,q}", "", 100, 50, 0., 1., "nph");
702  i_sector.second.m_correlationHistsY["ClusterProbXYQ"] = i_sector.second.bookCorrHistsY(
703  "ClusterProbXYQ", "cluster probability xyq", "prob_{cl,xyq}", "", 100, 50, 0., 1., "nph");
704  i_sector.second.m_correlationHistsY["LogClusterProb"] = i_sector.second.bookCorrHistsY(
705  "LogClusterProb", "cluster probability xy", "log(prob_{cl,xy})", "", 60, 30, logClusterProbMin, 0., "nph");
706  i_sector.second.m_correlationHistsY["IsOnEdge"] =
707  i_sector.second.bookCorrHistsY("IsOnEdge", "IsOnEdge", "isOnEdge", "", 2, 2, 0, 2, "nph");
708  i_sector.second.m_correlationHistsY["HasBadPixels"] =
709  i_sector.second.bookCorrHistsY("HasBadPixels", "HasBadPixels", "hasBadPixels", "", 2, 2, 0, 2, "nph");
710  i_sector.second.m_correlationHistsY["SpansTwoRoc"] =
711  i_sector.second.bookCorrHistsY("SpansTwoRoc", "SpansTwoRoc", "spansTwoRoc", "", 2, 2, 0, 2, "nph");
712  i_sector.second.m_correlationHistsY["QBin"] =
713  i_sector.second.bookCorrHistsY("QBin", "q bin", "q bin", "", 8, 8, 0, 8, "nph");
714  }
715 
716  else {
717  i_sector.second.m_correlationHistsX["ChargeStrip"] = i_sector.second.bookCorrHistsX(
718  "ChargeStrip", "cluster charge", "c_{cl}", "[APV counts]", 100, 50, 0., chargeStripMax, "nph");
719  i_sector.second.m_correlationHistsX["MaxStrip"] = i_sector.second.bookCorrHistsX(
720  "MaxStrip", "strip with max. charge", "n_{cl,max}", "[# strips]", 800, 800, -10., 790., "npht");
721  i_sector.second.m_correlationHistsX["MaxCharge"] = i_sector.second.bookCorrHistsX(
722  "MaxCharge", "charge of strip with max. charge", "c_{cl,max}", "[APV counts]", 300, 75, -10., 290., "nph");
723  i_sector.second.m_correlationHistsX["MaxIndex"] = i_sector.second.bookCorrHistsX(
724  "MaxIndex", "cluster-index of strip with max. charge", "i_{cl,max}", "[# strips]", 10, 10, 0., 10., "nph");
725  i_sector.second.m_correlationHistsX["ChargeOnEdges"] =
726  i_sector.second.bookCorrHistsX("ChargeOnEdges",
727  "fraction of charge on edge strips",
728  "(c_{st,L}+c_{st,R})/c_{cl}",
729  "",
730  60,
731  60,
732  -0.1,
733  1.1,
734  "nph");
735  i_sector.second.m_correlationHistsX["ChargeAsymmetry"] =
736  i_sector.second.bookCorrHistsX("ChargeAsymmetry",
737  "asymmetry of charge on edge strips",
738  "(c_{st,L}-c_{st,R})/c_{cl}",
739  "",
740  110,
741  55,
742  -1.1,
743  1.1,
744  "nph");
745  i_sector.second.m_correlationHistsX["ChargeLRplus"] =
746  i_sector.second.bookCorrHistsX("ChargeLRplus",
747  "fraction of charge not on maxStrip",
748  "(c_{cl,L}+c_{cl,R})/c_{cl}",
749  "",
750  60,
751  60,
752  -0.1,
753  1.1,
754  "nph");
755  i_sector.second.m_correlationHistsX["ChargeLRminus"] =
756  i_sector.second.bookCorrHistsX("ChargeLRminus",
757  "asymmetry of charge L and R of maxStrip",
758  "(c_{cl,L}-c_{cl,R})/c_{cl}",
759  "",
760  110,
761  55,
762  -1.1,
763  1.1,
764  "nph");
765  i_sector.second.m_correlationHistsX["SOverN"] =
766  i_sector.second.bookCorrHistsX("SOverN", "signal over noise", "s/N", "", 100, 50, 0, sOverNMax, "nph");
767  i_sector.second.m_correlationHistsX["WidthProj"] = i_sector.second.bookCorrHistsX(
768  "WidthProj", "projected width", "w_{p}", "[# strips]", 200, 20, 0., widthMax, "nph");
769  i_sector.second.m_correlationHistsX["WidthDiff"] = i_sector.second.bookCorrHistsX("WidthDiff",
770  "width difference",
771  "w_{p} - w_{cl}",
772  "[# strips]",
773  200,
774  20,
775  -widthMax / 2.,
776  widthMax / 2.,
777  "nph");
778 
779  i_sector.second.WidthVsWidthProjected = secDir.make<TH2F>("h2_widthVsWidthProj",
780  "w_{cl} vs. w_{p};w_{p} [# strips];w_{cl} [# strips]",
781  static_cast<int>(widthMax),
782  0,
783  widthMax,
784  static_cast<int>(widthMax),
785  0,
786  widthMax);
787  i_sector.second.PWidthVsWidthProjected =
788  secDir.make<TProfile>("p_widthVsWidthProj",
789  "w_{cl} vs. w_{p};w_{p} [# strips];w_{cl} [# strips]",
790  static_cast<int>(widthMax),
791  0,
792  widthMax);
793 
794  i_sector.second.WidthDiffVsMaxStrip =
795  secDir.make<TH2F>("h2_widthDiffVsMaxStrip",
796  "(w_{p} - w_{cl}) vs. n_{cl,max};n_{cl,max};w_{p} - w_{cl} [# strips]",
797  800,
798  -10.,
799  790.,
800  static_cast<int>(widthMax),
801  -widthMax / 2.,
802  widthMax / 2.);
803  i_sector.second.PWidthDiffVsMaxStrip =
804  secDir.make<TProfile>("p_widthDiffVsMaxStrip",
805  "(w_{p} - w_{cl}) vs. n_{cl,max};n_{cl,max};w_{p} - w_{cl} [# strips]",
806  800,
807  -10.,
808  790.);
809 
810  i_sector.second.WidthDiffVsSigmaXHit =
811  secDir.make<TH2F>("h2_widthDiffVsSigmaXHit",
812  "(w_{p} - w_{cl}) vs. #sigma_{hit,x};#sigma_{hit,x} [cm];w_{p} - w_{cl} [# strips]",
813  100,
814  0.,
815  sigmaXMax,
816  100,
817  -10.,
818  10.);
819  i_sector.second.PWidthDiffVsSigmaXHit =
820  secDir.make<TProfile>("p_widthDiffVsSigmaXHit",
821  "(w_{p} - w_{cl}) vs. #sigma_{hit,x};#sigma_{hit,x} [cm];w_{p} - w_{cl} [# strips]",
822  100,
823  0.,
824  sigmaXMax);
825 
826  i_sector.second.WidthVsPhiSensX =
827  secDir.make<TH2F>("h2_widthVsPhiSensX",
828  "w_{cl} vs. #phi_{module,x};#phi_{module,x} [ ^{o}];w_{cl} [# strips]",
829  93,
830  -93,
831  93,
832  static_cast<int>(widthMax),
833  0,
834  widthMax);
835  i_sector.second.PWidthVsPhiSensX = secDir.make<TProfile>(
836  "p_widthVsPhiSensX", "w_{cl} vs. #phi_{module,x};#phi_{module,x} [ ^{o}];w_{cl} [# strips]", 93, -93, 93);
837  }
838 
839  // Hit Parameters (transform errors and residuals from [cm] in [mum])
840  i_sector.second.m_correlationHistsX["SigmaXHit"] = i_sector.second.bookCorrHistsX(
841  "SigmaXHit", "hit error", "#sigma_{hit,x}", "[#mum]", 105, 20, sigmaXMin * 10000., sigmaXMax * 10000., "np");
842  i_sector.second.m_correlationHistsX["SigmaXTrk"] = i_sector.second.bookCorrHistsX(
843  "SigmaXTrk", "track error", "#sigma_{trk,x}", "[#mum]", 105, 20, sigmaXMin * 10000., sigmaXMax * 10000., "np");
844  i_sector.second.m_correlationHistsX["SigmaX"] = i_sector.second.bookCorrHistsX(
845  "SigmaX", "residual error", "#sigma_{r,x}", "[#mum]", 105, 20, sigmaXMin * 10000., sigmaXMax * 10000., "np");
846  i_sector.second.m_correlationHistsX["PhiSens"] = i_sector.second.bookCorrHistsX(
847  "PhiSens", "track angle on sensor", "#phi_{module}", "[ ^{o}]", 96, 48, -3, 93, "nphtr");
848  i_sector.second.m_correlationHistsX["PhiSensX"] = i_sector.second.bookCorrHistsX(
849  "PhiSensX", "track angle on sensor", "#phi_{module,x}", "[ ^{o}]", 186, 93, -phiSensXMax, phiSensXMax, "nphtr");
850  i_sector.second.m_correlationHistsX["PhiSensY"] = i_sector.second.bookCorrHistsX(
851  "PhiSensY", "track angle on sensor", "#phi_{module,y}", "[ ^{o}]", 186, 93, -93, 93, "nphtr");
852 
853  i_sector.second.XHit = secDir.make<TH1F>("h_XHit", " hit measurement x_{hit};x_{hit} [cm];# hits", 100, -20, 20);
854  i_sector.second.XTrk = secDir.make<TH1F>("h_XTrk", "track prediction x_{trk};x_{trk} [cm];# hits", 100, -20, 20);
855  i_sector.second.SigmaX2 =
856  secDir.make<TH1F>("h_SigmaX2",
857  "squared residual error #sigma_{r,x}^{2};#sigma_{r,x}^{2} [#mum^{2}];# hits",
858  105,
859  sigmaXMin * 10000.,
860  sigmaX2Max * 10000. * 10000.); //no mistake !
861  i_sector.second.ResX = secDir.make<TH1F>(
862  "h_ResX", "residual r_{x};x_{trk}-x_{hit} [#mum];# hits", 100, -resXAbsMax * 10000., resXAbsMax * 10000.);
863  i_sector.second.NorResX =
864  secDir.make<TH1F>("h_NorResX",
865  "normalized residual r_{x}/#sigma_{r,x};(x_{trk}-x_{hit})/#sigma_{r,x};# hits",
866  100,
867  -norResXAbsMax,
868  norResXAbsMax);
869  i_sector.second.ProbX = secDir.make<TH1F>(
870  "h_ProbX", "residual probability;prob(r_{x}^{2}/#sigma_{r,x}^{2},1);# hits", 60, probXMin, probXMax);
871 
872  i_sector.second.PhiSensXVsBarycentreX =
873  secDir.make<TH2F>("h2_phiSensXVsBarycentreX",
874  "#phi_{module,x} vs. b_{cl,x};b_{cl,x} [# channels];#phi_{module,x} [ ^{o}]",
875  200,
876  -10.,
877  790.,
878  93,
879  -93,
880  93);
881  i_sector.second.PPhiSensXVsBarycentreX =
882  secDir.make<TProfile>("p_phiSensXVsBarycentreX",
883  "#phi_{module,x} vs. b_{cl,x};b_{cl,x} [# channels];#phi_{module,x} [ ^{o}]",
884  200,
885  -10.,
886  790.);
887 
888  if (pixelSector) {
889  i_sector.second.m_correlationHistsY["SigmaYHit"] = i_sector.second.bookCorrHistsY(
890  "SigmaYHit", "hit error", "#sigma_{hit,y}", "[#mum]", 105, 20, sigmaXMin * 10000., sigmaXMax * 10000., "np");
891  i_sector.second.m_correlationHistsY["SigmaYTrk"] = i_sector.second.bookCorrHistsY(
892  "SigmaYTrk", "track error", "#sigma_{trk,y}", "[#mum]", 105, 20, sigmaXMin * 10000., sigmaXMax * 10000., "np");
893  i_sector.second.m_correlationHistsY["SigmaY"] = i_sector.second.bookCorrHistsY(
894  "SigmaY", "residual error", "#sigma_{r,y}", "[#mum]", 105, 20, sigmaXMin * 10000., sigmaXMax * 10000., "np");
895  i_sector.second.m_correlationHistsY["PhiSens"] = i_sector.second.bookCorrHistsY(
896  "PhiSens", "track angle on sensor", "#phi_{module}", "[ ^{o}]", 96, 48, -3, 93, "nphtr");
897  i_sector.second.m_correlationHistsY["PhiSensX"] = i_sector.second.bookCorrHistsY("PhiSensX",
898  "track angle on sensor",
899  "#phi_{module,x}",
900  "[ ^{o}]",
901  186,
902  93,
903  -phiSensXMax,
904  phiSensXMax,
905  "nphtr");
906  i_sector.second.m_correlationHistsY["PhiSensY"] = i_sector.second.bookCorrHistsY(
907  "PhiSensY", "track angle on sensor", "#phi_{module,y}", "[ ^{o}]", 186, 93, -93, 93, "nphtr");
908 
909  i_sector.second.YHit = secDir.make<TH1F>("h_YHit", " hit measurement y_{hit};y_{hit} [cm];# hits", 100, -20, 20);
910  i_sector.second.YTrk = secDir.make<TH1F>("h_YTrk", "track prediction y_{trk};y_{trk} [cm];# hits", 100, -20, 20);
911  i_sector.second.SigmaY2 =
912  secDir.make<TH1F>("h_SigmaY2",
913  "squared residual error #sigma_{r,y}^{2};#sigma_{r,y}^{2} [#mum^{2}];# hits",
914  105,
915  sigmaXMin * 10000.,
916  sigmaX2Max * 10000. * 10000.); //no mistake !
917  i_sector.second.ResY = secDir.make<TH1F>(
918  "h_ResY", "residual r_{y};y_{trk}-y_{hit} [#mum];# hits", 100, -resXAbsMax * 10000., resXAbsMax * 10000.);
919  i_sector.second.NorResY =
920  secDir.make<TH1F>("h_NorResY",
921  "normalized residual r_{y}/#sigma_{r,y};(y_{trk}-y_{hit})/#sigma_{r,y};# hits",
922  100,
923  -norResXAbsMax,
924  norResXAbsMax);
925  i_sector.second.ProbY = secDir.make<TH1F>(
926  "h_ProbY", "residual probability;prob(r_{y}^{2}/#sigma_{r,y}^{2},1);# hits", 60, probXMin, probXMax);
927 
928  i_sector.second.PhiSensYVsBarycentreY =
929  secDir.make<TH2F>("h2_phiSensYVsBarycentreY",
930  "#phi_{module,y} vs. b_{cl,y};b_{cl,y} [# channels];#phi_{module,y} [ ^{o}]",
931  200,
932  -10.,
933  790.,
934  93,
935  -93,
936  93);
937  i_sector.second.PPhiSensYVsBarycentreY =
938  secDir.make<TProfile>("p_phiSensYVsBarycentreY",
939  "#phi_{module,y} vs. b_{cl,y};b_{cl,y} [# channels];#phi_{module,y} [ ^{o}]",
940  200,
941  -10.,
942  790.);
943  }
944 
945  // Track Parameters
946  i_sector.second.m_correlationHistsX["HitsValid"] =
947  i_sector.second.bookCorrHistsX("HitsValid", "# hits", "[valid]", 50, 0, 50, "npt");
948  i_sector.second.m_correlationHistsX["HitsInvalid"] =
949  i_sector.second.bookCorrHistsX("HitsInvalid", "# hits", "[invalid]", 20, 0, 20, "npt");
950  i_sector.second.m_correlationHistsX["Hits2D"] =
951  i_sector.second.bookCorrHistsX("Hits2D", "# hits", "[2D]", 20, 0, 20, "npt");
952  i_sector.second.m_correlationHistsX["LayersMissed"] =
953  i_sector.second.bookCorrHistsX("LayersMissed", "# layers", "[missed]", 10, 0, 10, "npt");
954  i_sector.second.m_correlationHistsX["HitsPixel"] =
955  i_sector.second.bookCorrHistsX("HitsPixel", "# hits", "[pixel]", 10, 0, 10, "npt");
956  i_sector.second.m_correlationHistsX["HitsStrip"] =
957  i_sector.second.bookCorrHistsX("HitsStrip", "# hits", "[strip]", 40, 0, 40, "npt");
958  i_sector.second.m_correlationHistsX["HitsGood"] =
959  i_sector.second.bookCorrHistsX("HitsGood", "# hits", "[good]", 50, 0, 50, "npt");
960  i_sector.second.m_correlationHistsX["NorChi2"] =
961  i_sector.second.bookCorrHistsX("NorChi2", "#chi^{2}/f", "", 50, 0, norChi2Max, "npr");
962  i_sector.second.m_correlationHistsX["Theta"] =
963  i_sector.second.bookCorrHistsX("Theta", "#theta", "[ ^{o}]", 40, -10, 190, "npt");
964  i_sector.second.m_correlationHistsX["Phi"] =
965  i_sector.second.bookCorrHistsX("Phi", "#phi", "[ ^{o}]", 76, -190, 190, "npt");
966  i_sector.second.m_correlationHistsX["D0Beamspot"] =
967  i_sector.second.bookCorrHistsX("D0Beamspot", "d_{0, BS}", "[cm]", 40, -d0Max, d0Max, "npt");
968  i_sector.second.m_correlationHistsX["Dz"] =
969  i_sector.second.bookCorrHistsX("Dz", "d_{z}", "[cm]", 40, -dzMax, dzMax, "npt");
970  i_sector.second.m_correlationHistsX["Pt"] =
971  i_sector.second.bookCorrHistsX("Pt", "p_{t}", "[GeV]", 50, 0, pMax, "npt");
972  i_sector.second.m_correlationHistsX["P"] = i_sector.second.bookCorrHistsX("P", "|p|", "[GeV]", 50, 0, pMax, "npt");
973  i_sector.second.m_correlationHistsX["InvP"] =
974  i_sector.second.bookCorrHistsX("InvP", "1/|p|", "[GeV^{-1}]", 25, 0, invPMax, "t");
975  i_sector.second.m_correlationHistsX["MeanAngle"] =
976  i_sector.second.bookCorrHistsX("MeanAngle", "<#phi_{module}>", "[ ^{o}]", 25, -5, 95, "npt");
977  //i_sector.second.m_correlationHistsX[""] = i_sector.second.bookCorrHistsX("","","",,,,"nphtr");
978 
979  if (pixelSector) {
980  i_sector.second.m_correlationHistsY["HitsValid"] =
981  i_sector.second.bookCorrHistsY("HitsValid", "# hits", "[valid]", 50, 0, 50, "npt");
982  i_sector.second.m_correlationHistsY["HitsInvalid"] =
983  i_sector.second.bookCorrHistsY("HitsInvalid", "# hits", "[invalid]", 20, 0, 20, "npt");
984  i_sector.second.m_correlationHistsY["Hits2D"] =
985  i_sector.second.bookCorrHistsY("Hits2D", "# hits", "[2D]", 20, 0, 20, "npt");
986  i_sector.second.m_correlationHistsY["LayersMissed"] =
987  i_sector.second.bookCorrHistsY("LayersMissed", "# layers", "[missed]", 10, 0, 10, "npt");
988  i_sector.second.m_correlationHistsY["HitsPixel"] =
989  i_sector.second.bookCorrHistsY("HitsPixel", "# hits", "[pixel]", 10, 0, 10, "npt");
990  i_sector.second.m_correlationHistsY["HitsStrip"] =
991  i_sector.second.bookCorrHistsY("HitsStrip", "# hits", "[strip]", 40, 0, 40, "npt");
992  i_sector.second.m_correlationHistsY["HitsGood"] =
993  i_sector.second.bookCorrHistsY("HitsGood", "# hits", "[good]", 50, 0, 50, "npt");
994  i_sector.second.m_correlationHistsY["NorChi2"] =
995  i_sector.second.bookCorrHistsY("NorChi2", "#chi^{2}/f", "", 50, 0, norChi2Max, "npr");
996  i_sector.second.m_correlationHistsY["Theta"] =
997  i_sector.second.bookCorrHistsY("Theta", "#theta", "[ ^{o}]", 40, -10, 190, "npt");
998  i_sector.second.m_correlationHistsY["Phi"] =
999  i_sector.second.bookCorrHistsY("Phi", "#phi", "[ ^{o}]", 76, -190, 190, "npt");
1000  i_sector.second.m_correlationHistsY["D0Beamspot"] =
1001  i_sector.second.bookCorrHistsY("D0Beamspot", "d_{0, BS}", "[cm]", 40, -d0Max, d0Max, "npt");
1002  i_sector.second.m_correlationHistsY["Dz"] =
1003  i_sector.second.bookCorrHistsY("Dz", "d_{z}", "[cm]", 40, -dzMax, dzMax, "npt");
1004  i_sector.second.m_correlationHistsY["Pt"] =
1005  i_sector.second.bookCorrHistsY("Pt", "p_{t}", "[GeV]", 50, 0, pMax, "npt");
1006  i_sector.second.m_correlationHistsY["P"] =
1007  i_sector.second.bookCorrHistsY("P", "|p|", "[GeV]", 50, 0, pMax, "npt");
1008  i_sector.second.m_correlationHistsY["InvP"] =
1009  i_sector.second.bookCorrHistsY("InvP", "1/|p|", "[GeV^{-1}]", 25, 0, invPMax, "t");
1010  i_sector.second.m_correlationHistsY["MeanAngle"] =
1011  i_sector.second.bookCorrHistsY("MeanAngle", "<#phi_{module}>", "[ ^{o}]", 25, -5, 95, "npt");
1012  }
1013 
1014  // (transform errors and residuals from [cm] in [mum])
1015  for (auto const& i_errHists : v_errHists) {
1016  double xMin(0.01 * (i_errHists - 1)), xMax(0.01 * (i_errHists));
1017  std::stringstream sigmaXHit, sigmaXTrk, sigmaX;
1018  sigmaXHit << "h_sigmaXHit_" << i_errHists;
1019  sigmaXTrk << "h_sigmaXTrk_" << i_errHists;
1020  sigmaX << "h_sigmaX_" << i_errHists;
1021  i_sector.second.m_sigmaX["sigmaXHit"].push_back(
1022  secDir.make<TH1F>(sigmaXHit.str().c_str(),
1023  "hit error #sigma_{hit,x};#sigma_{hit,x} [#mum];# hits",
1024  100,
1025  xMin * 10000.,
1026  xMax * 10000.));
1027  i_sector.second.m_sigmaX["sigmaXTrk"].push_back(
1028  secDir.make<TH1F>(sigmaXTrk.str().c_str(),
1029  "track error #sigma_{trk,x};#sigma_{trk,x} [#mum];# hits",
1030  100,
1031  xMin * 10000.,
1032  xMax * 10000.));
1033  i_sector.second.m_sigmaX["sigmaX"].push_back(
1034  secDir.make<TH1F>(sigmaX.str().c_str(),
1035  "residual error #sigma_{r,x};#sigma_{r,x} [#mum];# hits",
1036  100,
1037  xMin * 10000.,
1038  xMax * 10000.));
1039  if (pixelSector) {
1040  std::stringstream sigmaYHit, sigmaYTrk, sigmaY;
1041  sigmaYHit << "h_sigmaYHit_" << i_errHists;
1042  sigmaYTrk << "h_sigmaYTrk_" << i_errHists;
1043  sigmaY << "h_sigmaY_" << i_errHists;
1044  i_sector.second.m_sigmaY["sigmaYHit"].push_back(
1045  secDir.make<TH1F>(sigmaYHit.str().c_str(),
1046  "hit error #sigma_{hit,y};#sigma_{hit,y} [#mum];# hits",
1047  100,
1048  xMin * 10000.,
1049  xMax * 10000.));
1050  i_sector.second.m_sigmaY["sigmaYTrk"].push_back(
1051  secDir.make<TH1F>(sigmaYTrk.str().c_str(),
1052  "track error #sigma_{trk,y};#sigma_{trk,y} [#mum];# hits",
1053  100,
1054  xMin * 10000.,
1055  xMax * 10000.));
1056  i_sector.second.m_sigmaY["sigmaY"].push_back(
1057  secDir.make<TH1F>(sigmaY.str().c_str(),
1058  "residual error #sigma_{r,y};#sigma_{r,y} [#mum];# hits",
1059  100,
1060  xMin * 10000.,
1061  xMax * 10000.));
1062  }
1063  }
1064  }
1065 }
1066 
1068  std::vector<unsigned int> v_errHists(parameterSet_.getParameter<std::vector<unsigned int> >("vErrHists"));
1069  for (std::vector<unsigned int>::iterator i_errHists = v_errHists.begin(); i_errHists != v_errHists.end();
1070  ++i_errHists) {
1071  for (std::vector<unsigned int>::iterator i_errHists2 = i_errHists; i_errHists2 != v_errHists.end();) {
1072  ++i_errHists2;
1073  if (*i_errHists == *i_errHists2) {
1074  edm::LogError("BookSectorHists") << "Value of vErrHists in config exists twice: " << *i_errHists
1075  << "\n... delete one of both";
1076  v_errHists.erase(i_errHists2);
1077  }
1078  }
1079  }
1080 
1081  for (auto& i_sector : m_tkSector_) {
1083  if (!fileService) {
1084  throw edm::Exception(edm::errors::Configuration, "TFileService is not registered in cfg file");
1085  }
1086 
1087  std::stringstream sector;
1088  sector << "Sector_" << i_sector.first;
1089  TFileDirectory secDir = fileService->mkdir(sector.str());
1090 
1091  // Dummy histo containing the sector name as title
1092  i_sector.second.Name = secDir.make<TH1F>("z_name", i_sector.second.name.c_str(), 1, 0, 1);
1093 
1094  // Do not book histos for empty sectors
1095  if (i_sector.second.v_rawId.empty()) {
1096  continue;
1097  }
1098 
1099  // Distributions in each interval (stay in [cm], to have all calculations in [cm])
1100  if (m_resErrBins_
1101  .empty()) { // default if no selection taken into account: calculate APE with one bin with residual error 0-100um
1102  m_resErrBins_[1].first = 0.;
1103  m_resErrBins_[1].second = 0.01;
1104  }
1105  for (auto const& i_errBins : m_resErrBins_) {
1106  std::stringstream interval;
1107  interval << "Interval_" << i_errBins.first;
1108  TFileDirectory intDir = secDir.mkdir(interval.str());
1109  i_sector.second.m_binnedHists[i_errBins.first]["sigmaX"] =
1110  intDir.make<TH1F>("h_sigmaX", "residual resolution #sigma_{x};#sigma_{x} [cm];# hits", 100, 0., 0.01);
1111  i_sector.second.m_binnedHists[i_errBins.first]["norResX"] = intDir.make<TH1F>(
1112  "h_norResX", "normalized residual r_{x}/#sigma_{r,x};(x_{trk}-x_{hit})/#sigma_{r,x};# hits", 100, -10, 10);
1113  if (i_sector.second.isPixel) {
1114  i_sector.second.m_binnedHists[i_errBins.first]["sigmaY"] =
1115  intDir.make<TH1F>("h_sigmaY", "residual resolution #sigma_{y};#sigma_{y} [cm];# hits", 100, 0., 0.01);
1116  i_sector.second.m_binnedHists[i_errBins.first]["norResY"] = intDir.make<TH1F>(
1117  "h_norResY", "normalized residual r_{y}/#sigma_{r,y};(y_{trk}-y_{hit})/#sigma_{r,y};# hits", 100, -10, 10);
1118  }
1119  }
1120 
1121  TFileDirectory resDir = secDir.mkdir("Results");
1122 
1123  // TTree containing rawIds of all modules in sector
1124  unsigned int rawId(0);
1125  i_sector.second.RawId = resDir.make<TTree>("rawIdTree", "Tree containing rawIds of all modules in sector");
1126  i_sector.second.RawId->Branch("RawId", &rawId, "RawId/i");
1127  for (auto const& i_rawId : i_sector.second.v_rawId) {
1128  rawId = i_rawId;
1129  i_sector.second.RawId->Fill();
1130  }
1131 
1132  // Result plots (one hist per sector containing one bin per interval)
1133  // (transform errors and residuals from [cm] in [mum])
1134  std::vector<double> v_binX(parameterSet_.getParameter<std::vector<double> >("residualErrorBinning"));
1135  for (auto& i_binX : v_binX) {
1136  i_binX *= 10000.;
1137  }
1138  i_sector.second.EntriesX =
1139  resDir.make<TH1F>("h_entriesX", "# hits used;#sigma_{x} [#mum];# hits", v_binX.size() - 1, &(v_binX[0]));
1140  if (i_sector.second.isPixel) {
1141  i_sector.second.EntriesY =
1142  resDir.make<TH1F>("h_entriesY", "# hits used;#sigma_{y} [#mum];# hits", v_binX.size() - 1, &(v_binX[0]));
1143  }
1144 
1145  // In fact these are un-needed Analyzer plots, but I want to have them always for every sector visible
1146  // (transform errors and residuals from [cm] in [mum])
1147  i_sector.second.ResX = resDir.make<TH1F>(
1148  "h_ResX", "residual r_{x};x_{trk}-x_{hit} [#mum];# hits", 100, -0.03 * 10000., 0.03 * 10000.);
1149  i_sector.second.NorResX = resDir.make<TH1F>(
1150  "h_NorResX", "normalized residual r_{x}/#sigma_{r,x};(x_{trk}-x_{hit})/#sigma_{r,x};# hits", 100, -5., 5.);
1151  if (i_sector.second.isPixel) {
1152  i_sector.second.ResY = resDir.make<TH1F>(
1153  "h_ResY", "residual r_{y};y_{trk}-y_{hit} [#mum];# hits", 100, -0.03 * 10000., 0.03 * 10000.);
1154  i_sector.second.NorResY = resDir.make<TH1F>(
1155  "h_NorResY", "normalized residual r_{y}/#sigma_{r,y};(y_{trk}-y_{hit})/#sigma_{r,y};# hits", 100, -5., 5.);
1156  }
1157  }
1158 }
1159 
1160 // -----------------------------------------------------------------------------------------------------------
1161 
1163  bool zoomHists(parameterSet_.getParameter<bool>("zoomHists"));
1164 
1165  int trackSizeBins = zoomHists ? 6 : 201;
1166  double trackSizeMax = trackSizeBins - 1;
1167 
1168  double chi2Max = zoomHists ? 100. : 2000.;
1169  double norChi2Max = zoomHists ? 5. : 1000.;
1170  double d0max = zoomHists ? 0.02 : 40.; // cosmics: 100.|100.
1171  double dzmax = zoomHists ? 15. : 100.; // cosmics: 200.|600.
1172  double pMax = zoomHists ? 200. : 2000.;
1173 
1175  TFileDirectory evtDir = fileService->mkdir("EventVariables");
1177  evtDir.make<TH1F>("h_trackSize", "# tracks [all];# tracks;# events", trackSizeBins, -1, trackSizeMax);
1179  evtDir.make<TH1F>("h_trackSizeGood", "# tracks [good];# tracks;# events", trackSizeBins, -1, trackSizeMax);
1180  TFileDirectory trkDir = fileService->mkdir("TrackVariables");
1181  tkDetector_.HitsSize = trkDir.make<TH1F>("h_hitsSize", "# hits;# hits;# tracks", 51, -1, 50);
1182  tkDetector_.HitsValid = trkDir.make<TH1F>("h_hitsValid", "# hits [valid];# hits [valid];# tracks", 51, -1, 50);
1184  trkDir.make<TH1F>("h_hitsInvalid", "# hits [invalid];# hits [invalid];# tracks", 21, -1, 20);
1185  tkDetector_.Hits2D = trkDir.make<TH1F>("h_hits2D", "# hits [2D];# hits [2D];# tracks", 21, -1, 20);
1187  trkDir.make<TH1F>("h_layersMissed", "# layers [missed];# layers [missed];# tracks", 11, -1, 10);
1188  tkDetector_.HitsPixel = trkDir.make<TH1F>("h_hitsPixel", "# hits [pixel];# hits [pixel];# tracks", 11, -1, 10);
1189  tkDetector_.HitsStrip = trkDir.make<TH1F>("h_hitsStrip", "# hits [strip];# hits [strip];# tracks", 41, -1, 40);
1190  tkDetector_.Charge = trkDir.make<TH1F>("h_charge", "charge q;q [e];# tracks", 5, -2, 3);
1191  tkDetector_.Chi2 = trkDir.make<TH1F>("h_chi2", " #chi^{2};#chi^{2};# tracks", 100, 0, chi2Max);
1192  tkDetector_.Ndof = trkDir.make<TH1F>("h_ndof", "# degrees of freedom f;f;# tracks", 101, -1, 100);
1193  tkDetector_.NorChi2 = trkDir.make<TH1F>("h_norChi2", "normalized #chi^{2};#chi^{2}/f;# tracks", 200, 0, norChi2Max);
1194  tkDetector_.Prob = trkDir.make<TH1F>("h_prob", " #chi^{2} probability;prob(#chi^{2},f);# tracks", 50, 0, 1);
1195  tkDetector_.Eta = trkDir.make<TH1F>("h_eta", "pseudorapidity #eta;#eta;# tracks", 100, -5, 5);
1196  tkDetector_.EtaErr = trkDir.make<TH1F>("h_etaErr", "Error of #eta;#sigma(#eta);# tracks", 100, 0, 0.001);
1198  trkDir.make<TH1F>("h_etaSig", "Significance of #eta;#eta/#sigma(#eta);# tracks", 100, -20000, 20000);
1199  tkDetector_.Theta = trkDir.make<TH1F>("h_theta", "polar angle #theta;#theta [ ^{o}];# tracks", 100, -10, 190);
1200  tkDetector_.Phi = trkDir.make<TH1F>("h_phi", "azimuth angle #phi;#phi [ ^{o}];# tracks", 190, -190, 190);
1201  tkDetector_.PhiErr = trkDir.make<TH1F>("h_phiErr", "Error of #phi;#sigma(#phi) [ ^{o}];# tracks", 100, 0, 0.04);
1203  trkDir.make<TH1F>("h_phiSig", "Significance of #phi;#phi/#sigma(#phi) [ ^{o}];# tracks", 100, -50000, 50000);
1204  tkDetector_.D0Beamspot = trkDir.make<TH1F>(
1205  "h_d0Beamspot", "Closest approach d_{0} wrt. beamspot;d_{0, BS} [cm];# tracks", 200, -d0max, d0max);
1207  trkDir.make<TH1F>("h_d0BeamspotErr", "Error of d_{0, BS};#sigma(d_{0, BS}) [cm];# tracks", 200, 0, 0.01);
1208  tkDetector_.D0BeamspotSig = trkDir.make<TH1F>(
1209  "h_d0BeamspotSig", "Significance of d_{0, BS};d_{0, BS}/#sigma(d_{0, BS});# tracks", 100, -5, 5);
1210  tkDetector_.Dz = trkDir.make<TH1F>("h_dz", "Closest approach d_{z};d_{z} [cm];# tracks", 200, -dzmax, dzmax);
1211  tkDetector_.DzErr = trkDir.make<TH1F>("h_dzErr", "Error of d_{z};#sigma(d_{z}) [cm];# tracks", 200, 0, 0.01);
1212  tkDetector_.DzSig =
1213  trkDir.make<TH1F>("h_dzSig", "Significance of d_{z};d_{z}/#sigma(d_{z});# tracks", 100, -10000, 10000);
1214  tkDetector_.Pt = trkDir.make<TH1F>("h_pt", "transverse momentum p_{t};p_{t} [GeV];# tracks", 100, 0, pMax);
1215  tkDetector_.PtErr = trkDir.make<TH1F>("h_ptErr", "Error of p_{t};#sigma(p_{t}) [GeV];# tracks", 100, 0, 1.6);
1216  tkDetector_.PtSig = trkDir.make<TH1F>("h_ptSig", "Significance of p_{t};p_{t}/#sigma(p_{t});# tracks", 100, 0, 200);
1217  tkDetector_.P = trkDir.make<TH1F>("h_p", "momentum magnitude |p|;|p| [GeV];# tracks", 100, 0, pMax);
1218  tkDetector_.MeanAngle = trkDir.make<TH1F>(
1219  "h_meanAngle", "mean angle on module <#phi_{module}>;<#phi_{module}> [ ^{o}];# tracks", 100, -5, 95);
1220  tkDetector_.HitsGood = trkDir.make<TH1F>("h_hitsGood", "# hits [good];# hits [good];# tracks", 51, -1, 50);
1221 
1222  tkDetector_.MeanAngleVsHits = trkDir.make<TH2F>(
1223  "h2_meanAngleVsHits", "<#phi_{module}> vs. # hits;# hits;<#phi_{module}> [ ^{o}]", 51, -1, 50, 50, -5, 95);
1225  trkDir.make<TH2F>("h2_hitsGoodVsHitsValid",
1226  "# hits [good] vs. # hits [valid];# hits [valid];# hits [good]",
1227  51,
1228  -1,
1229  50,
1230  51,
1231  -1,
1232  50);
1234  trkDir.make<TH2F>("h2_hitsPixelVsEta", "# hits [pixel] vs. #eta;#eta;# hits [pixel]", 60, -3, 3, 11, -1, 10);
1235  tkDetector_.HitsPixelVsTheta = trkDir.make<TH2F>(
1236  "h2_hitsPixelVsTheta", "# hits [pixel] vs. #theta;#theta;# hits [pixel]", 100, -10, 190, 11, -1, 10);
1238  trkDir.make<TH2F>("h2_hitsStripVsEta", "# hits [strip] vs. #eta;#eta;# hits [strip]", 60, -3, 3, 31, -1, 40);
1239  tkDetector_.HitsStripVsTheta = trkDir.make<TH2F>(
1240  "h2_hitsStripVsTheta", "# hits [strip] vs. #theta;#theta;# hits [strip]", 100, -10, 190, 31, -1, 40);
1241  tkDetector_.PtVsEta = trkDir.make<TH2F>("h2_ptVsEta", "p_{t} vs. #eta;#eta;p_{t} [GeV]", 60, -3, 3, 100, 0, pMax);
1243  trkDir.make<TH2F>("h2_ptVsTheta", "p_{t} vs. #theta;#theta;p_{t} [GeV]", 100, -10, 190, 100, 0, pMax);
1244 
1245  tkDetector_.PMeanAngleVsHits = trkDir.make<TProfile>(
1246  "p_meanAngleVsHits", "<#phi_{module}> vs. # hits;# hits;<#phi_{module}> [ ^{o}]", 51, -1, 50);
1247  tkDetector_.PHitsGoodVsHitsValid = trkDir.make<TProfile>(
1248  "p_hitsGoodVsHitsValid", "# hits [good] vs. # hits [valid];# hits [valid];# hits [good]", 51, -1, 50);
1250  trkDir.make<TProfile>("p_hitsPixelVsEta", "# hits [pixel] vs. #eta;#eta;# hits [pixel]", 60, -3, 3);
1252  trkDir.make<TProfile>("p_hitsPixelVsTheta", "# hits [pixel] vs. #theta;#theta;# hits [pixel]", 100, -10, 190);
1254  trkDir.make<TProfile>("p_hitsStripVsEta", "# hits [strip] vs. #eta;#eta;# hits [strip]", 60, -3, 3);
1256  trkDir.make<TProfile>("p_hitsStripVsTheta", "# hits [strip] vs. #theta;#theta;# hits [strip]", 100, -10, 190);
1257  tkDetector_.PPtVsEta = trkDir.make<TProfile>("p_ptVsEta", "p_{t} vs. #eta;#eta;p_{t} [GeV]", 60, -3, 3);
1258  tkDetector_.PPtVsTheta = trkDir.make<TProfile>("p_ptVsTheta", "p_{t} vs. #theta;#theta;p_{t} [GeV]", 100, -10, 190);
1259 }
1260 
1261 // -----------------------------------------------------------------------------------------------------------
1262 
1264  const Trajectory& traj,
1265  const reco::BeamSpot& beamSpot) {
1266  const math::XYZPoint beamPoint(beamSpot.x0(), beamSpot.y0(), beamSpot.z0());
1267  double d0BeamspotErr =
1268  std::sqrt(track.d0Error() * track.d0Error() + 0.5 * beamSpot.BeamWidthX() * beamSpot.BeamWidthX() +
1269  0.5 * beamSpot.BeamWidthY() * beamSpot.BeamWidthY());
1270 
1271  const static TrajectoryStateCombiner tsoscomb;
1272 
1273  const reco::HitPattern& hitPattern(track.hitPattern());
1274 
1276 
1277  trkParams.hitsSize = track.recHitsSize();
1278  trkParams.hitsValid = track.found(); // invalid is every hit from every single module that expects a hit
1279  trkParams.hitsInvalid = trkParams.hitsSize - trkParams.hitsValid;
1280  trkParams.layersMissed =
1281  track.lost(); // lost hit means, that a crossed layer doesn't contain a hit (can be more than one invalid hit)
1282  trkParams.hitsPixel = hitPattern.numberOfValidPixelHits();
1283  trkParams.hitsStrip = hitPattern.numberOfValidStripHits();
1284  trkParams.charge = track.charge();
1285  trkParams.chi2 = track.chi2();
1286  trkParams.ndof = track.ndof();
1287  trkParams.norChi2 = trkParams.chi2 / trkParams.ndof;
1288  trkParams.prob = TMath::Prob(trkParams.chi2, trkParams.ndof);
1289  trkParams.eta = track.eta();
1290  trkParams.etaErr = track.etaError();
1291  trkParams.theta = track.theta();
1292  trkParams.phi = track.phi();
1293  trkParams.phiErr = track.phiError();
1294  trkParams.d0 = track.d0();
1295  trkParams.d0Beamspot = -1. * track.dxy(beamPoint);
1296  trkParams.d0BeamspotErr = d0BeamspotErr;
1297  trkParams.dz = track.dz();
1298  trkParams.dzErr = track.dzError();
1299  trkParams.dzBeamspot = track.dz(beamPoint);
1300  trkParams.p = track.p();
1301  trkParams.pt = track.pt();
1302  trkParams.ptErr = track.ptError();
1303 
1304  const std::vector<TrajectoryMeasurement>& v_meas = traj.measurements();
1305 
1306  int count2D(0);
1307  float meanPhiSensToNorm(0.F);
1308  for (auto const& i_meas : v_meas) {
1309  const TrajectoryMeasurement& meas = i_meas;
1310  const TransientTrackingRecHit& hit = *meas.recHit();
1311  const TrackingRecHit& recHit = *hit.hit();
1312  if (this->isHit2D(recHit))
1313  ++count2D;
1314 
1316  const align::LocalVector mom(tsos.localDirection());
1317  meanPhiSensToNorm += atan(fabs(sqrt(mom.x() * mom.x() + mom.y() * mom.y()) / mom.z()));
1318  }
1319  meanPhiSensToNorm *= (1. / static_cast<float>(trkParams.hitsSize));
1320 
1321  trkParams.hits2D = count2D;
1322  trkParams.meanPhiSensToNorm = meanPhiSensToNorm;
1323 
1324  if (parameterSet_.getParameter<bool>("applyTrackCuts")) {
1325  trackCut_ = false;
1326  if (trkParams.hitsStrip < 11 || trkParams.hits2D < 2 || trkParams.hitsPixel < 2 || //trkParams.hitsInvalid>2 ||
1327  trkParams.hitsStrip > 35 || trkParams.hitsPixel > 7 || trkParams.norChi2 > 5. || trkParams.pt < 25. ||
1328  trkParams.pt > 150. || std::abs(trkParams.d0Beamspot) > 0.02 || std::abs(trkParams.dz) > 15.)
1329  trackCut_ = true;
1330  } else {
1331  trackCut_ = false;
1332  }
1333 
1334  return trkParams;
1335 }
1336 
1338  const edm::EventSetup& iSetup) {
1340 
1341  const static TrajectoryStateCombiner tsoscomb;
1342 
1343  const TrajectoryMeasurement& meas = i_meas;
1344  const TransientTrackingRecHit& hit = *meas.recHit();
1345  const TrackingRecHit& recHit = *hit.hit();
1346  const TrajectoryStateOnSurface& tsos = tsoscomb(meas.forwardPredictedState(), meas.backwardPredictedState());
1347 
1348  const DetId& detId(hit.geographicalId());
1349  const DetId::Detector& detector = detId.det();
1350  if (detector != DetId::Tracker) {
1351  hitParams.hitState = TrackStruct::notInTracker;
1352  return hitParams;
1353  }
1354  const uint32_t rawId(detId.rawId());
1355 
1356  for (auto& i_sector : m_tkSector_) {
1357  for (auto const& i_rawId : i_sector.second.v_rawId) {
1358  if (rawId == i_rawId) {
1359  hitParams.v_sector.push_back(i_sector.first);
1360  break;
1361  }
1362  }
1363  }
1364 
1365  const align::LocalVector& mom(tsos.localDirection());
1366  int xMomentum(0), yMomentum(0), zMomentum(0);
1367  xMomentum = mom.x() > 0. ? 1 : -1;
1368  yMomentum = mom.y() > 0. ? 1 : -1;
1369  zMomentum = mom.z() > 0. ? 1 : -1;
1370  float phiSensX =
1371  std::atan(std::fabs(mom.x() / mom.z())) *
1372  static_cast<float>(
1373  m_tkTreeVar_[rawId].vDirection); // check for orientation of E- and B- Field (thoughts for barrel)
1374  float phiSensY = std::atan(std::fabs(mom.y() / mom.z())) * static_cast<float>(m_tkTreeVar_[rawId].vDirection);
1375  hitParams.phiSens = std::atan(std::fabs(std::sqrt(mom.x() * mom.x() + mom.y() * mom.y()) / mom.z()));
1376  hitParams.phiSensX = (xMomentum == zMomentum ? phiSensX : -phiSensX);
1377  hitParams.phiSensY = (yMomentum == zMomentum ? phiSensY : -phiSensY);
1378 
1379  if (!hit.isValid()) {
1380  hitParams.hitState = TrackStruct::invalid;
1381  return hitParams;
1382  }
1383 
1384  // Get local positions and errors of hit and track
1385  const LocalPoint& lPHit = hit.localPosition();
1386  const LocalPoint& lPTrk = tsos.localPosition();
1387 
1388  // use APE also for the hit error, while APE is automatically included in tsos error
1389  //
1390  // no need to add APE to hitError anymore by Ajay 27 Oct 2014
1391  const LocalError& errHitApe = hit.localPositionError(); // now sum of CPE+APE as said by MARCO?
1392  LocalError errorWithoutAPE;
1393 
1394  bool Pixel(false);
1395  bool Strip(false);
1396 
1397  if (m_tkTreeVar_[rawId].subdetId == PixelSubdetector::PixelBarrel ||
1398  m_tkTreeVar_[rawId].subdetId == PixelSubdetector::PixelEndcap) {
1399  Pixel = true;
1400  } else if (m_tkTreeVar_[rawId].subdetId == StripSubdetector::TIB ||
1401  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TOB ||
1402  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TID ||
1403  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TEC) {
1404  Strip = true;
1405  } else {
1406  edm::LogWarning("FillHitVariables") << "cant identify wether hit is from pixel or strip";
1407  hitParams.hitState = TrackStruct::invalid;
1408  return hitParams;
1409  }
1410 
1411  if (!hit.detUnit()) {
1412  hitParams.hitState = TrackStruct::invalid;
1413  return hitParams;
1414  } // is it a single physical module?
1415  const GeomDetUnit& detUnit = *hit.detUnit();
1416 
1417  if (Pixel) {
1418  if (!dynamic_cast<const PixelTopology*>(&detUnit.type().topology())) {
1419  hitParams.hitState = TrackStruct::invalid;
1420  return hitParams;
1421  }
1422  const PixelGeomDetUnit* pixelDet = (const PixelGeomDetUnit*)(&detUnit);
1423  const LocalError& lape = pixelDet->localAlignmentError();
1424  if (lape.valid()) {
1425  errorWithoutAPE = LocalError(errHitApe.xx() - lape.xx(), errHitApe.xy() - lape.xy(), errHitApe.yy() - lape.yy());
1426  }
1427  }
1428  if (Strip) {
1429  if (!dynamic_cast<const StripTopology*>(&detUnit.type().topology())) {
1430  hitParams.hitState = TrackStruct::invalid;
1431  return hitParams;
1432  }
1433  const StripGeomDetUnit* stripDet = (const StripGeomDetUnit*)(&detUnit);
1434  const LocalError& lape = stripDet->localAlignmentError();
1435  if (lape.valid()) {
1436  errorWithoutAPE = LocalError(errHitApe.xx() - lape.xx(), errHitApe.xy() - lape.xy(), errHitApe.yy() - lape.yy());
1437  }
1438  }
1439 
1440  const LocalError& errHitWoApe = errorWithoutAPE;
1441  const LocalError& errTrk = tsos.localError().positionError();
1442 
1443  const StatePositionAndError2 positionAndError2Hit = this->positionAndError2(lPHit, errHitApe, hit);
1444  const StatePositionAndError2 positionAndError2HitWoApe = this->positionAndError2(lPHit, errHitWoApe, hit);
1445  edm::LogInfo("CalculateAPE") << "errHitWoApe " << errHitWoApe << "errHitApe " << errHitApe;
1446 
1447  const StatePositionAndError2 positionAndError2Trk = this->positionAndError2(lPTrk, errTrk, hit);
1448 
1449  const TrackStruct::HitState& stateHit(positionAndError2Hit.first);
1450  const TrackStruct::HitState& stateHitWoApe(positionAndError2HitWoApe.first);
1451  const TrackStruct::HitState& stateTrk(positionAndError2Trk.first);
1452 
1453  if (stateHit == TrackStruct::invalid || stateHitWoApe == TrackStruct::invalid || stateTrk == TrackStruct::invalid) {
1454  hitParams.hitState = TrackStruct::invalid;
1455  return hitParams;
1456  } else if (stateHit == TrackStruct::negativeError || stateHitWoApe == TrackStruct::negativeError ||
1457  stateTrk == TrackStruct::negativeError) {
1458  ++counter1;
1459  // Do not print error message by default
1460  //std::stringstream ss_error;
1461  //ss_error<<"Upper values belong to: ";
1462  //if(stateHit==TrackStruct::negativeError)ss_error<<"Hit without APE, ";
1463  //if(stateHitWoApe==TrackStruct::negativeError)ss_error<<"Hit with APE, ";
1464  //if(stateTrk==TrackStruct::negativeError)ss_error<<"Track,";
1465  //edm::LogError("Negative error Value")<<"@SUB=ApeEstimator::fillHitVariables"<<ss_error.str();
1467  return hitParams;
1468  }
1469 
1470  // Calculate residuals
1471 
1472  const float xHit = positionAndError2Hit.second.posX;
1473  const float xTrk = positionAndError2Trk.second.posX;
1474  const float yHit = positionAndError2Hit.second.posY;
1475  const float yTrk = positionAndError2Trk.second.posY;
1476 
1477  const float errXHit2(positionAndError2Hit.second.errX2);
1478  const float errXHitWoApe2(positionAndError2HitWoApe.second.errX2);
1479  const float errXTrk2(positionAndError2Trk.second.errX2);
1480  const float errYHit2(positionAndError2Hit.second.errY2);
1481  const float errYHitWoApe2(positionAndError2HitWoApe.second.errY2);
1482  const float errYTrk2(positionAndError2Trk.second.errY2);
1483 
1484  const float errXHit = std::sqrt(positionAndError2Hit.second.errX2);
1485  const float errXHitWoApe = std::sqrt(positionAndError2HitWoApe.second.errX2);
1486  const float errXTrk = std::sqrt(positionAndError2Trk.second.errX2);
1487  const float errYHit = std::sqrt(positionAndError2Hit.second.errY2);
1488  const float errYHitWoApe = std::sqrt(positionAndError2HitWoApe.second.errY2);
1489  const float errYTrk = std::sqrt(positionAndError2Trk.second.errY2);
1490 
1491  const float resX = xTrk - xHit;
1492  const float resY = yTrk - yHit;
1493 
1494  const float errX = std::sqrt(errXHit2 + errXTrk2);
1495  const float errXWoApe2 = errXHitWoApe2 + errXTrk2;
1496  const float errXWoApe = std::sqrt(errXWoApe2);
1497  const float errY = std::sqrt(errYHit2 + errYTrk2);
1498  const float errYWoApe2 = errYHitWoApe2 + errYTrk2;
1499  const float errYWoApe = std::sqrt(errYWoApe2);
1500 
1501  const float norResX = resX / errX;
1502  const float norResY = resY / errY;
1503 
1504  // Take global orientation into account for residuals (sign is not important for errors)
1505 
1506  float resXprime(999.F), resYprime(999.F), norResXprime(999.F), norResYprime(999.F);
1507  if (m_tkTreeVar_[rawId].uDirection == 1) {
1508  resXprime = resX;
1509  norResXprime = norResX;
1510  } else if (m_tkTreeVar_[rawId].uDirection == -1) {
1511  resXprime = -resX;
1512  norResXprime = -norResX;
1513  } else {
1514  edm::LogError("FillHitVariables") << "Incorrect value of uDirection, which gives global module orientation";
1515  hitParams.hitState = TrackStruct::invalid;
1516  return hitParams;
1517  }
1518  if (m_tkTreeVar_[rawId].vDirection == 1) {
1519  resYprime = resY;
1520  norResYprime = norResY;
1521  } else if (m_tkTreeVar_[rawId].vDirection == -1) {
1522  resYprime = -resY;
1523  norResYprime = -norResY;
1524  } else {
1525  edm::LogError("FillHitVariables") << "Incorrect value of vDirection, which gives global module orientation";
1526  hitParams.hitState = TrackStruct::invalid;
1527  return hitParams;
1528  }
1529 
1530  hitParams.xHit = xHit;
1531  hitParams.xTrk = xTrk;
1532 
1533  hitParams.errXHit = errXHit;
1534  hitParams.errXHitWoApe = errXHitWoApe;
1535  hitParams.errXTrk = errXTrk;
1536 
1537  hitParams.errX2 = errX * errX;
1538  hitParams.errX = errX;
1539  hitParams.errXWoApe = errXWoApe;
1540 
1541  hitParams.resX = resXprime;
1542  hitParams.norResX = norResXprime;
1543 
1544  const float norResX2(norResXprime * norResXprime);
1545  hitParams.probX = TMath::Prob(norResX2, 1);
1546 
1547  hitParams.yHit = yHit;
1548  hitParams.yTrk = yTrk;
1549 
1550  hitParams.errYHit = errYHit;
1551  hitParams.errYHitWoApe = errYHitWoApe;
1552  hitParams.errYTrk = errYTrk;
1553 
1554  hitParams.errY2 = errY * errY;
1555  hitParams.errY = errY;
1556  hitParams.errYWoApe = errYWoApe;
1557 
1558  hitParams.resY = resYprime;
1559  hitParams.norResY = norResYprime;
1560 
1561  const float norResY2(norResYprime * norResYprime);
1562  hitParams.probY = TMath::Prob(norResY2, 1);
1563 
1564  // Cluster parameters
1565 
1566  if (m_tkTreeVar_[rawId].subdetId == PixelSubdetector::PixelBarrel ||
1567  m_tkTreeVar_[rawId].subdetId == PixelSubdetector::PixelEndcap) {
1568  const SiPixelRecHit& pixelHit = dynamic_cast<const SiPixelRecHit&>(recHit);
1569  const SiPixelCluster& pixelCluster = *pixelHit.cluster();
1570 
1571  hitParams.chargePixel = pixelCluster.charge();
1572  hitParams.widthX = pixelCluster.sizeX();
1573  hitParams.baryStripX = pixelCluster.x();
1574  hitParams.widthY = pixelCluster.sizeY();
1575  hitParams.baryStripY = pixelCluster.y();
1576 
1577  hitParams.clusterProbabilityXY = pixelHit.clusterProbability(0);
1578  hitParams.clusterProbabilityQ = pixelHit.clusterProbability(2);
1579  hitParams.clusterProbabilityXYQ = pixelHit.clusterProbability(1);
1580  hitParams.logClusterProbability = std::log10(hitParams.clusterProbabilityXY);
1581 
1582  hitParams.isOnEdge = pixelHit.isOnEdge();
1583  hitParams.hasBadPixels = pixelHit.hasBadPixels();
1584  hitParams.spansTwoRoc = pixelHit.spansTwoROCs();
1585  hitParams.qBin = pixelHit.qBin();
1586 
1587  hitParams.isPixelHit = true;
1588  } else if (m_tkTreeVar_[rawId].subdetId == StripSubdetector::TIB ||
1589  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TOB ||
1590  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TID ||
1591  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TEC) {
1592  if (!(dynamic_cast<const SiStripRecHit2D*>(&recHit) || dynamic_cast<const SiStripRecHit1D*>(&recHit))) {
1593  edm::LogError("FillHitVariables")
1594  << "RecHit in Strip is 'Matched' or 'Projected', but here all should be monohits per module";
1595  hitParams.hitState = TrackStruct::invalid;
1596  return hitParams;
1597  }
1598  const SiStripCluster* clusterPtr(nullptr);
1599  if (m_tkTreeVar_[rawId].subdetId == StripSubdetector::TIB ||
1600  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TOB) {
1601  if (dynamic_cast<const SiStripRecHit1D*>(&recHit)) {
1602  const SiStripRecHit1D& stripHit = dynamic_cast<const SiStripRecHit1D&>(recHit);
1603  clusterPtr = &(*stripHit.cluster());
1604  } else if (dynamic_cast<const SiStripRecHit2D*>(&recHit)) {
1605  edm::LogWarning("FillHitVariables") << "Data has TIB/TOB hits as SiStripRecHit2D and not 1D. Probably data is "
1606  "processed with CMSSW<34X. Nevertheless everything should work fine";
1607  const SiStripRecHit2D& stripHit = dynamic_cast<const SiStripRecHit2D&>(recHit);
1608  clusterPtr = &(*stripHit.cluster());
1609  }
1610  } else if (m_tkTreeVar_[rawId].subdetId == StripSubdetector::TID ||
1611  m_tkTreeVar_[rawId].subdetId == StripSubdetector::TEC) {
1612  const SiStripRecHit2D& stripHit = dynamic_cast<const SiStripRecHit2D&>(recHit);
1613  clusterPtr = &(*stripHit.cluster());
1614  }
1615  if (!clusterPtr) {
1616  edm::LogError("FillHitVariables") << "Pointer to cluster not valid!!! This should never happen...";
1617  hitParams.hitState = TrackStruct::invalid;
1618  return hitParams;
1619  }
1620  const SiStripCluster& stripCluster(*clusterPtr);
1621 
1622  const SiStripClusterInfo clusterInfo = SiStripClusterInfo(stripCluster, iSetup, rawId, std::string(""));
1623 
1624  const std::vector<uint8_t>::const_iterator stripChargeL(clusterInfo.stripCharges().begin());
1625  const std::vector<uint8_t>::const_iterator stripChargeR(--(clusterInfo.stripCharges().end()));
1626  const std::pair<uint16_t, uint16_t> stripChargeLR = std::make_pair(*stripChargeL, *stripChargeR);
1627 
1628  hitParams.chargeStrip = clusterInfo.charge();
1629  hitParams.widthX = clusterInfo.width();
1630  hitParams.baryStripX = clusterInfo.baryStrip() + 1.;
1631  hitParams.isModuleUsable = clusterInfo.IsModuleUsable();
1632  hitParams.maxStrip = clusterInfo.maxStrip() + 1;
1633  hitParams.maxStripInv = m_tkTreeVar_[rawId].nStrips - hitParams.maxStrip + 1;
1634  hitParams.maxCharge = clusterInfo.maxCharge();
1635  hitParams.maxIndex = clusterInfo.maxIndex();
1636  hitParams.chargeOnEdges =
1637  static_cast<float>(stripChargeLR.first + stripChargeLR.second) / static_cast<float>(hitParams.chargeStrip);
1638  hitParams.chargeAsymmetry = static_cast<float>(stripChargeLR.first - stripChargeLR.second) /
1639  static_cast<float>(stripChargeLR.first + stripChargeLR.second);
1640  hitParams.chargeLRplus = static_cast<float>(clusterInfo.chargeLR().first + clusterInfo.chargeLR().second) /
1641  static_cast<float>(hitParams.chargeStrip);
1642  hitParams.chargeLRminus = static_cast<float>(clusterInfo.chargeLR().first - clusterInfo.chargeLR().second) /
1643  static_cast<float>(hitParams.chargeStrip);
1644  hitParams.sOverN = clusterInfo.signalOverNoise();
1645 
1646  // Calculate projection length corrected by drift
1647  if (!hit.detUnit()) {
1648  hitParams.hitState = TrackStruct::invalid;
1649  return hitParams;
1650  } // is it a single physical module?
1651  const GeomDetUnit& detUnit = *hit.detUnit();
1652  if (!dynamic_cast<const StripTopology*>(&detUnit.type().topology())) {
1653  hitParams.hitState = TrackStruct::invalid;
1654  return hitParams;
1655  }
1656 
1657  edm::ESHandle<MagneticField> magFieldHandle;
1658  iSetup.get<IdealMagneticFieldRecord>().get(magFieldHandle);
1659 
1660  edm::ESHandle<SiStripLorentzAngle> lorentzAngleHandle;
1661  iSetup.get<SiStripLorentzAngleDepRcd>().get(lorentzAngleHandle); //MODIFIED BY LOIC QUERTENMONT
1662 
1663  const StripGeomDetUnit* stripDet = (const StripGeomDetUnit*)(&detUnit);
1664  const MagneticField* magField(magFieldHandle.product());
1665  LocalVector bField = (stripDet->surface()).toLocal(magField->inTesla(stripDet->surface().position()));
1666  const SiStripLorentzAngle* lorentzAngle(lorentzAngleHandle.product());
1667  float tanLorentzAnglePerTesla = lorentzAngle->getLorentzAngle(stripDet->geographicalId().rawId());
1668 
1669  float dirX = -tanLorentzAnglePerTesla * bField.y();
1670  float dirY = tanLorentzAnglePerTesla * bField.x();
1671  float dirZ = 1.; // E field always in z direction
1672  LocalVector driftDirection(dirX, dirY, dirZ);
1673 
1674  const Bounds& bounds = stripDet->specificSurface().bounds();
1675  float maxLength = std::sqrt(std::pow(bounds.length(), 2) + std::pow(bounds.width(), 2));
1676  float thickness = bounds.thickness();
1677 
1678  const StripTopology& topol = dynamic_cast<const StripTopology&>(detUnit.type().topology());
1679  LocalVector momentumDir(tsos.localDirection());
1680  LocalPoint momentumPos(tsos.localPosition());
1681  LocalVector scaledMomentumDir(momentumDir);
1682  if (momentumDir.z() > 0.)
1683  scaledMomentumDir *= std::fabs(thickness / momentumDir.z());
1684  else if (momentumDir.z() < 0.)
1685  scaledMomentumDir *= -std::fabs(thickness / momentumDir.z());
1686  else
1687  scaledMomentumDir *= maxLength / momentumDir.mag();
1688 
1689  float projEdge1 = topol.measurementPosition(momentumPos - 0.5 * scaledMomentumDir).x();
1690  if (projEdge1 < 0.)
1691  projEdge1 = 0.;
1692  else if (projEdge1 > m_tkTreeVar_[rawId].nStrips)
1693  projEdge1 = m_tkTreeVar_[rawId].nStrips;
1694  float projEdge2 = topol.measurementPosition(momentumPos + 0.5 * scaledMomentumDir).x();
1695  if (projEdge2 < 0.)
1696  projEdge1 = 0.;
1697  else if (projEdge2 > m_tkTreeVar_[rawId].nStrips)
1698  projEdge1 = m_tkTreeVar_[rawId].nStrips;
1699 
1700  float coveredStrips = std::fabs(projEdge2 - projEdge1);
1701 
1702  hitParams.projWidth = coveredStrips;
1703 
1704  } else {
1705  edm::LogError("FillHitVariables") << "Incorrect subdetector ID, hit not associated to tracker";
1706  hitParams.hitState = TrackStruct::notInTracker;
1707  return hitParams;
1708  }
1709 
1710  if (!hitParams.isModuleUsable) {
1711  hitParams.hitState = TrackStruct::invalid;
1712  return hitParams;
1713  }
1714 
1715  if (hitParams.v_sector.empty()) {
1717  return hitParams;
1718  }
1719 
1720  return hitParams;
1721  //}
1722 }
1723 
1725  const LocalError& localError,
1726  const TransientTrackingRecHit& hit) {
1728 
1729  const DetId& detId(hit.geographicalId());
1730  const uint32_t& rawId(detId.rawId());
1731  const unsigned int& subdetId(m_tkTreeVar_[rawId].subdetId);
1732 
1733  if (localError.xx() < 0. || localError.yy() < 0.) {
1734  // Do not print error message by default
1735  //edm::LogError("Negative error Value")<<"@SUB=ApeEstimator::fillHitVariables"
1736  // <<"One of the squared error methods gives negative result\n"
1737  // <<"\tSubdetector\tlocalError.xx()\tlocalError.yy()\n"
1738  // <<"\t"<<subdetId<<"\t\t"<<localError.xx()<<"\t"<<localError.yy();
1739  vPE2.first = TrackStruct::negativeError;
1740  return vPE2;
1741  }
1742 
1743  if (subdetId == PixelSubdetector::PixelBarrel || subdetId == PixelSubdetector::PixelEndcap ||
1744  subdetId == StripSubdetector::TIB || subdetId == StripSubdetector::TOB) {
1745  // Cartesian coordinates
1746  vPE2 = std::make_pair(TrackStruct::ok, this->rectangularPositionAndError2(localPoint, localError));
1747  } else if (subdetId == StripSubdetector::TID || subdetId == StripSubdetector::TEC) {
1748  // Local x in radial coordinates
1749  if (!hit.detUnit())
1750  return vPE2; // is it a single physical module?
1751  const GeomDetUnit& detUnit = *hit.detUnit();
1752 
1753  if (!dynamic_cast<const RadialStripTopology*>(&detUnit.type().topology()))
1754  return vPE2;
1755  const RadialStripTopology& topol = dynamic_cast<const RadialStripTopology&>(detUnit.type().topology());
1756 
1757  MeasurementError measError = topol.measurementError(localPoint, localError);
1758  if (measError.uu() < 0. || measError.vv() < 0.) {
1759  // Do not print error message by default
1760  //edm::LogError("Negative error Value")<<"@SUB=ApeEstimator::fillHitVariables"
1761  // <<"One of the squared error methods gives negative result\n"
1762  // <<"\tmeasError.uu()\tmeasError.vv()\n"
1763  // <<"\t"<<measError.uu()<<"\t"<<measError.vv()
1764  // <<"\n\nOriginalValues:\n"
1765  // <<localPoint.x()<<" "<<localPoint.y()<<"\n"
1766  // <<localError.xx()<<" "<<localError.yy()<<"\n"
1767  // <<"Subdet: "<<subdetId;
1768  vPE2.first = TrackStruct::negativeError;
1769  return vPE2;
1770  }
1771  vPE2 = std::make_pair(TrackStruct::ok, this->radialPositionAndError2(localPoint, localError, topol));
1772  } else {
1773  edm::LogError("FillHitVariables") << "Incorrect subdetector ID, hit not associated to tracker";
1774  }
1775 
1776  return vPE2;
1777 }
1778 
1780  const float x(lP.x());
1781  const float y(lP.y());
1782  const float errX2(lE.xx());
1783  const float errY2(lE.yy());
1784 
1785  return PositionAndError2(x, y, errX2, errY2);
1786 }
1787 
1789  const LocalError& lE,
1790  const RadialStripTopology& topol) {
1791  MeasurementPoint measPos = topol.measurementPosition(lP);
1792  MeasurementError measErr = topol.measurementError(lP, lE);
1793 
1794  const float r_0 = topol.originToIntersection();
1795  const float stripLength = topol.localStripLength(lP);
1796  const float phi = topol.stripAngle(measPos.x());
1797 
1798  float x(-999.F);
1799  float y(-999.F);
1800  float errX2(-999.F);
1801  float errY2(-999.F);
1802 
1803  x = phi * r_0;
1804  // Radial y (not symmetric around 0; radial distance with minimum at middle strip at lower edge [0, yMax])
1805  const float l_0 = r_0 - topol.detHeight() / 2;
1806  const float cosPhi(std::cos(phi));
1807  y = measPos.y() * stripLength - 0.5 * stripLength + l_0 * (1. / cosPhi - 1.);
1808 
1809  const float angularWidth2(topol.angularWidth() * topol.angularWidth());
1810  const float errPhi2(measErr.uu() * angularWidth2);
1811 
1812  errX2 = errPhi2 * r_0 * r_0;
1813  // Radial y (not symmetric around 0, real radial distance from intersection point)
1814  const float cosPhi4(std::pow(cosPhi, 4)), sinPhi2(std::sin(phi) * std::sin(phi));
1815  const float helpSummand = l_0 * l_0 * (sinPhi2 / cosPhi4 * errPhi2);
1816  errY2 = measErr.vv() * stripLength * stripLength + helpSummand;
1817 
1818  return PositionAndError2(x, y, errX2, errY2);
1819 }
1820 
1821 // -----------------------------------------------------------------------------------------------------------
1822 
1824  this->setHitSelectionMapUInt("width");
1825  this->setHitSelectionMap("widthProj");
1826  this->setHitSelectionMap("widthDiff");
1827  this->setHitSelectionMap("charge");
1828  this->setHitSelectionMapUInt("edgeStrips");
1829  this->setHitSelectionMap("maxCharge");
1830  this->setHitSelectionMapUInt("maxIndex");
1831  this->setHitSelectionMap("chargeOnEdges");
1832  this->setHitSelectionMap("chargeAsymmetry");
1833  this->setHitSelectionMap("chargeLRplus");
1834  this->setHitSelectionMap("chargeLRminus");
1835  this->setHitSelectionMap("sOverN");
1836 
1837  this->setHitSelectionMap("chargePixel");
1838  this->setHitSelectionMapUInt("widthX");
1839  this->setHitSelectionMapUInt("widthY");
1840 
1841  this->setHitSelectionMap("baryStripX");
1842  this->setHitSelectionMap("baryStripY");
1843  this->setHitSelectionMap("clusterProbabilityXY");
1844  this->setHitSelectionMap("clusterProbabilityQ");
1845  this->setHitSelectionMap("clusterProbabilityXYQ");
1846  this->setHitSelectionMap("logClusterProbability");
1847  this->setHitSelectionMapUInt("isOnEdge");
1848  this->setHitSelectionMapUInt("hasBadPixels");
1849  this->setHitSelectionMapUInt("spansTwoRoc");
1850  this->setHitSelectionMapUInt("qBin");
1851 
1852  this->setHitSelectionMap("phiSens");
1853  this->setHitSelectionMap("phiSensX");
1854  this->setHitSelectionMap("phiSensY");
1855  this->setHitSelectionMap("resX");
1856  this->setHitSelectionMap("norResX");
1857  this->setHitSelectionMap("probX");
1858  this->setHitSelectionMap("errXHit");
1859  this->setHitSelectionMap("errXTrk");
1860  this->setHitSelectionMap("errX");
1861  this->setHitSelectionMap("errX2");
1862 
1863  this->setHitSelectionMap("resY");
1864  this->setHitSelectionMap("norResY");
1865  this->setHitSelectionMap("probY");
1866  this->setHitSelectionMap("errYHit");
1867  this->setHitSelectionMap("errYTrk");
1868  this->setHitSelectionMap("errY");
1869  this->setHitSelectionMap("errY2");
1870 
1871  edm::LogInfo("HitSelector") << "applying hit cuts ...";
1872  bool emptyMap(true);
1873  for (auto& i_hitSelection : m_hitSelection_) {
1874  if (!i_hitSelection.second.empty()) {
1875  int entry(1);
1876  double intervalBegin(999.);
1877  for (std::vector<double>::iterator i_hitInterval = i_hitSelection.second.begin();
1878  i_hitInterval != i_hitSelection.second.end();
1879  ++entry) {
1880  if (entry % 2 == 1) {
1881  intervalBegin = *i_hitInterval;
1882  ++i_hitInterval;
1883  } else {
1884  if (intervalBegin > *i_hitInterval) {
1885  edm::LogError("HitSelector") << "INVALID Interval selected for " << i_hitSelection.first << ":\t"
1886  << intervalBegin << " > " << (*i_hitInterval) << "\n ... delete Selection for "
1887  << i_hitSelection.first;
1888  i_hitSelection.second.clear();
1889  i_hitInterval = i_hitSelection.second.begin(); //emptyMap = true; i_hitSelection = m_hitSelection_.begin();
1890  } else {
1891  edm::LogInfo("HitSelector") << "Interval selected for " << i_hitSelection.first << ":\t" << intervalBegin
1892  << ", " << (*i_hitInterval);
1893  ++i_hitInterval;
1894  }
1895  }
1896  }
1897  if (!i_hitSelection.second.empty())
1898  emptyMap = false;
1899  }
1900  }
1901 
1902  bool emptyMapUInt(true);
1903  for (auto& i_hitSelection : m_hitSelectionUInt_) {
1904  if (!i_hitSelection.second.empty()) {
1905  int entry(1);
1906  unsigned int intervalBegin(999);
1907  for (std::vector<unsigned int>::iterator i_hitInterval = i_hitSelection.second.begin();
1908  i_hitInterval != i_hitSelection.second.end();
1909  ++entry) {
1910  if (entry % 2 == 1) {
1911  intervalBegin = *i_hitInterval;
1912  ++i_hitInterval;
1913  } else {
1914  if (intervalBegin > *i_hitInterval) {
1915  edm::LogError("HitSelector") << "INVALID Interval selected for " << i_hitSelection.first << ":\t"
1916  << intervalBegin << " > " << (*i_hitInterval) << "\n ... delete Selection for "
1917  << i_hitSelection.first;
1918  i_hitSelection.second.clear();
1919  i_hitInterval = i_hitSelection.second.begin(); //emptyMap = true; i_hitSelection = m_hitSelection_.begin();
1920  } else {
1921  edm::LogInfo("HitSelector") << "Interval selected for " << i_hitSelection.first << ":\t" << intervalBegin
1922  << ", " << (*i_hitInterval);
1923  ++i_hitInterval;
1924  }
1925  }
1926  }
1927  if (!i_hitSelection.second.empty())
1928  emptyMapUInt = false;
1929  }
1930  }
1931 
1932  if (emptyMap && emptyMapUInt) {
1933  m_hitSelection_.clear();
1934  m_hitSelectionUInt_.clear();
1935  edm::LogInfo("HitSelector") << "NO hit cuts applied";
1936  }
1937  return;
1938 }
1939 
1942  std::vector<double> v_cutVariable(parSet.getParameter<std::vector<double> >(cutVariable));
1943  if (v_cutVariable.size() % 2 == 1) {
1944  edm::LogError("HitSelector") << "Invalid Hit Selection for " << cutVariable
1945  << ": need even number of arguments (intervals)"
1946  << "\n ... delete Selection for " << cutVariable;
1947  v_cutVariable.clear();
1948  m_hitSelection_[cutVariable] = v_cutVariable;
1949  return;
1950  }
1951  m_hitSelection_[cutVariable] = v_cutVariable;
1952  return;
1953 }
1954 
1957  std::vector<unsigned int> v_cutVariable(parSet.getParameter<std::vector<unsigned int> >(cutVariable));
1958  if (v_cutVariable.size() % 2 == 1) {
1959  edm::LogError("HitSelector") << "Invalid Hit Selection for " << cutVariable
1960  << ": need even number of arguments (intervals)"
1961  << "\n ... delete Selection for " << cutVariable;
1962  v_cutVariable.clear();
1963  m_hitSelectionUInt_[cutVariable] = v_cutVariable;
1964  return;
1965  }
1966  m_hitSelectionUInt_[cutVariable] = v_cutVariable;
1967  return;
1968 }
1969 
1970 // -----------------------------------------------------------------------------------------------------------
1971 
1973  if (hitParams.hitState == TrackStruct::notInTracker)
1974  return false;
1975  if (hitParams.hitState == TrackStruct::invalid || hitParams.hitState == TrackStruct::negativeError)
1976  return false;
1977 
1978  bool isGoodHit(true);
1979  bool isGoodHitX(true);
1980  bool isGoodHitY(true);
1981 
1982  for (auto& i_hitSelection : m_hitSelection_) {
1983  const std::string& hitSelection(i_hitSelection.first);
1984  const std::vector<double>& v_hitSelection(i_hitSelection.second);
1985  if (v_hitSelection.empty())
1986  continue;
1987 
1988  // For pixel and strip sectors in common
1989  if (hitSelection == "phiSens") {
1990  if (!this->inDoubleInterval(v_hitSelection, hitParams.phiSens))
1991  isGoodHit = false;
1992  } else if (hitSelection == "phiSensX") {
1993  if (!this->inDoubleInterval(v_hitSelection, hitParams.phiSensX))
1994  isGoodHit = false;
1995  } else if (hitSelection == "phiSensY") {
1996  if (!this->inDoubleInterval(v_hitSelection, hitParams.phiSensY))
1997  isGoodHit = false;
1998  }
1999 
2000  else if (hitSelection == "resX") {
2001  if (!this->inDoubleInterval(v_hitSelection, hitParams.resX))
2002  isGoodHitX = false;
2003  } else if (hitSelection == "norResX") {
2004  if (!this->inDoubleInterval(v_hitSelection, hitParams.norResX))
2005  isGoodHitX = false;
2006  } else if (hitSelection == "probX") {
2007  if (!this->inDoubleInterval(v_hitSelection, hitParams.probX))
2008  isGoodHitX = false;
2009  } else if (hitSelection == "errXHit") {
2010  if (!this->inDoubleInterval(v_hitSelection, hitParams.errXHit))
2011  isGoodHitX = false;
2012  } else if (hitSelection == "errXTrk") {
2013  if (!this->inDoubleInterval(v_hitSelection, hitParams.errXTrk))
2014  isGoodHitX = false;
2015  } else if (hitSelection == "errX") {
2016  if (!this->inDoubleInterval(v_hitSelection, hitParams.errX))
2017  isGoodHitX = false;
2018  } else if (hitSelection == "errX2") {
2019  if (!this->inDoubleInterval(v_hitSelection, hitParams.errX2))
2020  isGoodHitX = false;
2021  }
2022 
2023  // For pixel only
2024  if (hitParams.isPixelHit) {
2025  if (hitSelection == "chargePixel") {
2026  if (!this->inDoubleInterval(v_hitSelection, hitParams.chargePixel))
2027  isGoodHit = false;
2028  } else if (hitSelection == "clusterProbabilityXY") {
2029  if (!this->inDoubleInterval(v_hitSelection, hitParams.clusterProbabilityXY))
2030  isGoodHit = false;
2031  } else if (hitSelection == "clusterProbabilityQ") {
2032  if (!this->inDoubleInterval(v_hitSelection, hitParams.clusterProbabilityQ))
2033  isGoodHit = false;
2034  } else if (hitSelection == "clusterProbabilityXYQ") {
2035  if (!this->inDoubleInterval(v_hitSelection, hitParams.clusterProbabilityXYQ))
2036  isGoodHit = false;
2037  } else if (hitSelection == "logClusterProbability") {
2038  if (!this->inDoubleInterval(v_hitSelection, hitParams.logClusterProbability))
2039  isGoodHit = false;
2040  }
2041 
2042  else if (hitSelection == "baryStripX") {
2043  if (!this->inDoubleInterval(v_hitSelection, hitParams.baryStripX))
2044  isGoodHitX = false;
2045  } else if (hitSelection == "baryStripY") {
2046  if (!this->inDoubleInterval(v_hitSelection, hitParams.baryStripY))
2047  isGoodHitY = false;
2048  }
2049 
2050  else if (hitSelection == "resY") {
2051  if (!this->inDoubleInterval(v_hitSelection, hitParams.resY))
2052  isGoodHitY = false;
2053  } else if (hitSelection == "norResY") {
2054  if (!this->inDoubleInterval(v_hitSelection, hitParams.norResY))
2055  isGoodHitY = false;
2056  } else if (hitSelection == "probY") {
2057  if (!this->inDoubleInterval(v_hitSelection, hitParams.probY))
2058  isGoodHitY = false;
2059  } else if (hitSelection == "errYHit") {
2060  if (!this->inDoubleInterval(v_hitSelection, hitParams.errYHit))
2061  isGoodHitY = false;
2062  } else if (hitSelection == "errYTrk") {
2063  if (!this->inDoubleInterval(v_hitSelection, hitParams.errYTrk))
2064  isGoodHitY = false;
2065  } else if (hitSelection == "errY") {
2066  if (!this->inDoubleInterval(v_hitSelection, hitParams.errY))
2067  isGoodHitY = false;
2068  } else if (hitSelection == "errY2") {
2069  if (!this->inDoubleInterval(v_hitSelection, hitParams.errY2))
2070  isGoodHitY = false;
2071  }
2072  } else { // For strip only
2073  if (hitSelection == "widthProj") {
2074  if (!this->inDoubleInterval(v_hitSelection, hitParams.projWidth))
2075  isGoodHit = false;
2076  } else if (hitSelection == "widthDiff") {
2077  if (!this->inDoubleInterval(v_hitSelection, hitParams.projWidth - static_cast<float>(hitParams.widthX)))
2078  isGoodHit = false;
2079  } else if (hitSelection == "charge") {
2080  if (!this->inDoubleInterval(v_hitSelection, hitParams.chargeStrip))
2081  isGoodHit = false;
2082  } else if (hitSelection == "maxCharge") {
2083  if (!this->inDoubleInterval(v_hitSelection, hitParams.maxCharge))
2084  isGoodHit = false;
2085  } else if (hitSelection == "chargeOnEdges") {
2086  if (!this->inDoubleInterval(v_hitSelection, hitParams.chargeOnEdges))
2087  isGoodHit = false;
2088  } else if (hitSelection == "chargeAsymmetry") {
2089  if (!this->inDoubleInterval(v_hitSelection, hitParams.chargeAsymmetry))
2090  isGoodHit = false;
2091  } else if (hitSelection == "chargeLRplus") {
2092  if (!this->inDoubleInterval(v_hitSelection, hitParams.chargeLRplus))
2093  isGoodHit = false;
2094  } else if (hitSelection == "chargeLRminus") {
2095  if (!this->inDoubleInterval(v_hitSelection, hitParams.chargeLRminus))
2096  isGoodHit = false;
2097  } else if (hitSelection == "sOverN") {
2098  if (!this->inDoubleInterval(v_hitSelection, hitParams.sOverN))
2099  isGoodHit = false;
2100  }
2101  }
2102  }
2103 
2104  for (auto& i_hitSelection : m_hitSelectionUInt_) {
2105  const std::string& hitSelection(i_hitSelection.first);
2106  const std::vector<unsigned int>& v_hitSelection(i_hitSelection.second);
2107  if (v_hitSelection.empty())
2108  continue;
2109 
2110  // For pixel and strip sectors in common
2111 
2112  // For pixel only
2113  if (hitParams.isPixelHit) {
2114  if (hitSelection == "isOnEdge") {
2115  if (!this->inUintInterval(v_hitSelection, hitParams.isOnEdge))
2116  isGoodHit = false;
2117  } else if (hitSelection == "hasBadPixels") {
2118  if (!this->inUintInterval(v_hitSelection, hitParams.hasBadPixels))
2119  isGoodHit = false;
2120  } else if (hitSelection == "spansTwoRoc") {
2121  if (!this->inUintInterval(v_hitSelection, hitParams.spansTwoRoc))
2122  isGoodHit = false;
2123  } else if (hitSelection == "qBin") {
2124  if (!this->inUintInterval(v_hitSelection, hitParams.qBin))
2125  isGoodHit = false;
2126  }
2127 
2128  else if (hitSelection == "widthX") {
2129  if (!this->inUintInterval(v_hitSelection, hitParams.widthX))
2130  isGoodHitX = false;
2131  } else if (hitSelection == "widthY") {
2132  if (!this->inUintInterval(v_hitSelection, hitParams.widthY))
2133  isGoodHitY = false;
2134  }
2135  } else { // For strip only
2136  if (hitSelection == "width") {
2137  if (!this->inUintInterval(v_hitSelection, hitParams.widthX))
2138  isGoodHit = false;
2139  } else if (hitSelection == "edgeStrips") {
2140  if (!this->inUintInterval(v_hitSelection, hitParams.maxStrip, hitParams.maxStripInv))
2141  isGoodHit = false;
2142  } else if (hitSelection == "maxIndex") {
2143  if (!this->inUintInterval(v_hitSelection, hitParams.maxIndex))
2144  isGoodHit = false;
2145  }
2146  }
2147  }
2148 
2149  if (hitParams.isPixelHit) {
2150  hitParams.goodXMeasurement = isGoodHit && isGoodHitX;
2151  hitParams.goodYMeasurement = isGoodHit && isGoodHitY;
2152  } else {
2153  hitParams.goodXMeasurement = isGoodHit && isGoodHitX;
2154  hitParams.goodYMeasurement = false;
2155  }
2156 
2157  if (!hitParams.goodXMeasurement && !hitParams.goodYMeasurement)
2158  return false;
2159  else
2160  return true;
2161 }
2162 
2163 bool ApeEstimator::inDoubleInterval(const std::vector<double>& v_hitSelection, const float variable) const {
2164  int entry(0);
2165  double intervalBegin(999.);
2166  bool isSelected(false);
2167  for (auto const& i_hitInterval : v_hitSelection) {
2168  ++entry;
2169  if (entry % 2 == 1)
2170  intervalBegin = i_hitInterval;
2171  else if (variable >= intervalBegin && variable < i_hitInterval)
2172  isSelected = true;
2173  }
2174  return isSelected;
2175 }
2176 
2177 bool ApeEstimator::inUintInterval(const std::vector<unsigned int>& v_hitSelection,
2178  const unsigned int variable,
2179  const unsigned int variable2) const {
2180  int entry(0);
2181  unsigned int intervalBegin(999);
2182  bool isSelected(false);
2183  for (auto i_hitInterval : v_hitSelection) {
2184  ++entry;
2185  if (entry % 2 == 1)
2186  intervalBegin = i_hitInterval;
2187  else if (variable >= intervalBegin && variable <= i_hitInterval) {
2188  if (variable2 == 999 || (variable2 >= intervalBegin && variable2 <= i_hitInterval))
2189  isSelected = true;
2190  }
2191  }
2192  return isSelected;
2193 }
2194 
2195 // -----------------------------------------------------------------------------------------------------------
2196 
2198  unsigned int goodHitsPerTrack(trackStruct.v_hitParams.size());
2199  tkDetector_.HitsGood->Fill(goodHitsPerTrack);
2200  tkDetector_.HitsGoodVsHitsValid->Fill(trackStruct.trkParams.hitsValid, goodHitsPerTrack);
2201  tkDetector_.PHitsGoodVsHitsValid->Fill(trackStruct.trkParams.hitsValid, goodHitsPerTrack);
2202 
2203  if (parameterSet_.getParameter<bool>("applyTrackCuts")) {
2204  // which tracks to take? need min. nr. of selected hits?
2205  if (goodHitsPerTrack < minGoodHitsPerTrack_)
2206  return;
2207  }
2208 
2209  tkDetector_.HitsSize->Fill(trackStruct.trkParams.hitsSize);
2210  tkDetector_.HitsValid->Fill(trackStruct.trkParams.hitsValid);
2211  tkDetector_.HitsInvalid->Fill(trackStruct.trkParams.hitsInvalid);
2212  tkDetector_.Hits2D->Fill(trackStruct.trkParams.hits2D);
2213  tkDetector_.LayersMissed->Fill(trackStruct.trkParams.layersMissed);
2214  tkDetector_.HitsPixel->Fill(trackStruct.trkParams.hitsPixel);
2215  tkDetector_.HitsStrip->Fill(trackStruct.trkParams.hitsStrip);
2216  tkDetector_.Charge->Fill(trackStruct.trkParams.charge);
2217  tkDetector_.Chi2->Fill(trackStruct.trkParams.chi2);
2218  tkDetector_.Ndof->Fill(trackStruct.trkParams.ndof);
2219  tkDetector_.NorChi2->Fill(trackStruct.trkParams.norChi2);
2220  tkDetector_.Prob->Fill(trackStruct.trkParams.prob);
2221  tkDetector_.Eta->Fill(trackStruct.trkParams.eta);
2222  tkDetector_.EtaErr->Fill(trackStruct.trkParams.etaErr);
2223  tkDetector_.EtaSig->Fill(trackStruct.trkParams.eta / trackStruct.trkParams.etaErr);
2224  tkDetector_.Theta->Fill(trackStruct.trkParams.theta * 180. / M_PI);
2225  tkDetector_.Phi->Fill(trackStruct.trkParams.phi * 180. / M_PI);
2226  tkDetector_.PhiErr->Fill(trackStruct.trkParams.phiErr * 180. / M_PI);
2227  tkDetector_.PhiSig->Fill(trackStruct.trkParams.phi / trackStruct.trkParams.phiErr);
2228  tkDetector_.D0Beamspot->Fill(trackStruct.trkParams.d0Beamspot);
2229  tkDetector_.D0BeamspotErr->Fill(trackStruct.trkParams.d0BeamspotErr);
2230  tkDetector_.D0BeamspotSig->Fill(trackStruct.trkParams.d0Beamspot / trackStruct.trkParams.d0BeamspotErr);
2231  tkDetector_.Dz->Fill(trackStruct.trkParams.dz);
2232  tkDetector_.DzErr->Fill(trackStruct.trkParams.dzErr);
2233  tkDetector_.DzSig->Fill(trackStruct.trkParams.dz / trackStruct.trkParams.dzErr);
2234  tkDetector_.P->Fill(trackStruct.trkParams.p);
2235  tkDetector_.Pt->Fill(trackStruct.trkParams.pt);
2236  tkDetector_.PtErr->Fill(trackStruct.trkParams.ptErr);
2237  tkDetector_.PtSig->Fill(trackStruct.trkParams.pt / trackStruct.trkParams.ptErr);
2238  tkDetector_.MeanAngle->Fill(trackStruct.trkParams.meanPhiSensToNorm * 180. / M_PI);
2239 
2240  tkDetector_.MeanAngleVsHits->Fill(trackStruct.trkParams.hitsSize,
2241  trackStruct.trkParams.meanPhiSensToNorm * 180. / M_PI);
2242  tkDetector_.HitsPixelVsEta->Fill(trackStruct.trkParams.eta, trackStruct.trkParams.hitsPixel);
2243  tkDetector_.HitsPixelVsTheta->Fill(trackStruct.trkParams.theta * 180. / M_PI, trackStruct.trkParams.hitsPixel);
2244  tkDetector_.HitsStripVsEta->Fill(trackStruct.trkParams.eta, trackStruct.trkParams.hitsStrip);
2245  tkDetector_.HitsStripVsTheta->Fill(trackStruct.trkParams.theta * 180. / M_PI, trackStruct.trkParams.hitsStrip);
2246  tkDetector_.PtVsEta->Fill(trackStruct.trkParams.eta, trackStruct.trkParams.pt);
2247  tkDetector_.PtVsTheta->Fill(trackStruct.trkParams.theta * 180. / M_PI, trackStruct.trkParams.pt);
2248 
2249  tkDetector_.PMeanAngleVsHits->Fill(trackStruct.trkParams.hitsSize,
2250  trackStruct.trkParams.meanPhiSensToNorm * 180. / M_PI);
2251  tkDetector_.PHitsPixelVsEta->Fill(trackStruct.trkParams.eta, trackStruct.trkParams.hitsPixel);
2252  tkDetector_.PHitsPixelVsTheta->Fill(trackStruct.trkParams.theta * 180. / M_PI, trackStruct.trkParams.hitsPixel);
2253  tkDetector_.PHitsStripVsEta->Fill(trackStruct.trkParams.eta, trackStruct.trkParams.hitsStrip);
2254  tkDetector_.PHitsStripVsTheta->Fill(trackStruct.trkParams.theta * 180. / M_PI, trackStruct.trkParams.hitsStrip);
2255  tkDetector_.PPtVsEta->Fill(trackStruct.trkParams.eta, trackStruct.trkParams.pt);
2256  tkDetector_.PPtVsTheta->Fill(trackStruct.trkParams.theta * 180. / M_PI, trackStruct.trkParams.pt);
2257 
2258  for (auto& i_hit : trackStruct.v_hitParams) {
2259  const TrackStruct::HitParameterStruct& hit(i_hit);
2260  //Put here from earlier method
2262  continue;
2263 
2264  for (auto& i_sector : m_tkSector_) {
2265  bool moduleInSector(false);
2266  for (auto const& i_hitSector : hit.v_sector) {
2267  if (i_sector.first == i_hitSector) {
2268  moduleInSector = true;
2269  break;
2270  }
2271  }
2272  if (!moduleInSector)
2273  continue;
2274  TrackerSectorStruct& sector(i_sector.second);
2275 
2276  if (hit.goodXMeasurement) {
2277  std::map<std::string, TrackerSectorStruct::CorrelationHists>& m_corrHists(sector.m_correlationHistsX);
2278 
2279  // Cluster and Hit Parameters
2280  this->fillHitHistsXForAnalyzerMode(hit, sector);
2281 
2282  // Track Parameters
2283  m_corrHists["HitsValid"].fillCorrHistsX(hit, trackStruct.trkParams.hitsValid);
2284  m_corrHists["HitsGood"].fillCorrHistsX(hit, goodHitsPerTrack);
2285  m_corrHists["HitsInvalid"].fillCorrHistsX(hit, trackStruct.trkParams.hitsInvalid);
2286  m_corrHists["Hits2D"].fillCorrHistsX(hit, trackStruct.trkParams.hits2D);
2287  m_corrHists["LayersMissed"].fillCorrHistsX(hit, trackStruct.trkParams.layersMissed);
2288  m_corrHists["HitsPixel"].fillCorrHistsX(hit, trackStruct.trkParams.hitsPixel);
2289  m_corrHists["HitsStrip"].fillCorrHistsX(hit, trackStruct.trkParams.hitsStrip);
2290  m_corrHists["NorChi2"].fillCorrHistsX(hit, trackStruct.trkParams.norChi2);
2291  m_corrHists["Theta"].fillCorrHistsX(hit, trackStruct.trkParams.theta * 180. / M_PI);
2292  m_corrHists["Phi"].fillCorrHistsX(hit, trackStruct.trkParams.phi * 180. / M_PI);
2293  m_corrHists["D0Beamspot"].fillCorrHistsX(hit, trackStruct.trkParams.d0Beamspot);
2294  m_corrHists["Dz"].fillCorrHistsX(hit, trackStruct.trkParams.dz);
2295  m_corrHists["Pt"].fillCorrHistsX(hit, trackStruct.trkParams.pt);
2296  m_corrHists["P"].fillCorrHistsX(hit, trackStruct.trkParams.p);
2297  m_corrHists["InvP"].fillCorrHistsX(hit, 1. / trackStruct.trkParams.p);
2298  m_corrHists["MeanAngle"].fillCorrHistsX(hit, trackStruct.trkParams.meanPhiSensToNorm * 180. / M_PI);
2299  }
2300 
2301  if (hit.goodYMeasurement) {
2302  std::map<std::string, TrackerSectorStruct::CorrelationHists>& m_corrHists(sector.m_correlationHistsY);
2303 
2304  // Cluster and Hit Parameters
2305  this->fillHitHistsYForAnalyzerMode(hit, sector);
2306 
2307  // Track Parameters
2308  m_corrHists["HitsValid"].fillCorrHistsY(hit, trackStruct.trkParams.hitsValid);
2309  m_corrHists["HitsGood"].fillCorrHistsY(hit, goodHitsPerTrack);
2310  m_corrHists["HitsInvalid"].fillCorrHistsY(hit, trackStruct.trkParams.hitsInvalid);
2311  m_corrHists["Hits2D"].fillCorrHistsY(hit, trackStruct.trkParams.hits2D);
2312  m_corrHists["LayersMissed"].fillCorrHistsY(hit, trackStruct.trkParams.layersMissed);
2313  m_corrHists["HitsPixel"].fillCorrHistsY(hit, trackStruct.trkParams.hitsPixel);
2314  m_corrHists["HitsStrip"].fillCorrHistsY(hit, trackStruct.trkParams.hitsStrip);
2315  m_corrHists["NorChi2"].fillCorrHistsY(hit, trackStruct.trkParams.norChi2);
2316  m_corrHists["Theta"].fillCorrHistsY(hit, trackStruct.trkParams.theta * 180. / M_PI);
2317  m_corrHists["Phi"].fillCorrHistsY(hit, trackStruct.trkParams.phi * 180. / M_PI);
2318  m_corrHists["D0Beamspot"].fillCorrHistsY(hit, trackStruct.trkParams.d0Beamspot);
2319  m_corrHists["Dz"].fillCorrHistsY(hit, trackStruct.trkParams.dz);
2320  m_corrHists["Pt"].fillCorrHistsY(hit, trackStruct.trkParams.pt);
2321  m_corrHists["P"].fillCorrHistsY(hit, trackStruct.trkParams.p);
2322  m_corrHists["InvP"].fillCorrHistsY(hit, 1. / trackStruct.trkParams.p);
2323  m_corrHists["MeanAngle"].fillCorrHistsY(hit, trackStruct.trkParams.meanPhiSensToNorm * 180. / M_PI);
2324  }
2325 
2326  // Special Histograms
2327  for (auto& i_sigmaX : sector.m_sigmaX) {
2328  for (auto& iHist : i_sigmaX.second) {
2329  if (i_sigmaX.first == "sigmaXHit")
2330  iHist->Fill(hit.errXHit * 10000.);
2331  else if (i_sigmaX.first == "sigmaXTrk")
2332  iHist->Fill(hit.errXTrk * 10000.);
2333  else if (i_sigmaX.first == "sigmaX")
2334  iHist->Fill(hit.errX * 10000.);
2335  }
2336  }
2337  for (auto& i_sigmaY : sector.m_sigmaY) {
2338  for (auto& iHist : i_sigmaY.second) {
2339  if (i_sigmaY.first == "sigmaYHit")
2340  iHist->Fill(hit.errYHit * 10000.);
2341  else if (i_sigmaY.first == "sigmaYTrk")
2342  iHist->Fill(hit.errYTrk * 10000.);
2343  else if (i_sigmaY.first == "sigmaY")
2344  iHist->Fill(hit.errY * 10000.);
2345  }
2346  }
2347  }
2348  }
2349 }
2350 
2352  TrackerSectorStruct& sector) {
2353  std::map<std::string, TrackerSectorStruct::CorrelationHists>& m_corrHists(sector.m_correlationHistsX);
2354 
2355  // Cluster Parameters
2356  m_corrHists["WidthX"].fillCorrHistsX(hit, hit.widthX);
2357  m_corrHists["BaryStripX"].fillCorrHistsX(hit, hit.baryStripX);
2358 
2359  if (hit.isPixelHit) {
2360  m_corrHists["ChargePixel"].fillCorrHistsX(hit, hit.chargePixel);
2361  m_corrHists["ClusterProbXY"].fillCorrHistsX(hit, hit.clusterProbabilityXY);
2362  m_corrHists["ClusterProbQ"].fillCorrHistsX(hit, hit.clusterProbabilityQ);
2363  m_corrHists["ClusterProbXYQ"].fillCorrHistsX(hit, hit.clusterProbabilityXYQ);
2364  m_corrHists["LogClusterProb"].fillCorrHistsX(hit, hit.logClusterProbability);
2365  m_corrHists["IsOnEdge"].fillCorrHistsX(hit, hit.isOnEdge);
2366  m_corrHists["HasBadPixels"].fillCorrHistsX(hit, hit.hasBadPixels);
2367  m_corrHists["SpansTwoRoc"].fillCorrHistsX(hit, hit.spansTwoRoc);
2368  m_corrHists["QBin"].fillCorrHistsX(hit, hit.qBin);
2369 
2370  } else {
2371  m_corrHists["ChargeStrip"].fillCorrHistsX(hit, hit.chargeStrip);
2372  m_corrHists["MaxStrip"].fillCorrHistsX(hit, hit.maxStrip);
2373  m_corrHists["MaxCharge"].fillCorrHistsX(hit, hit.maxCharge);
2374  m_corrHists["MaxIndex"].fillCorrHistsX(hit, hit.maxIndex);
2375  m_corrHists["ChargeOnEdges"].fillCorrHistsX(hit, hit.chargeOnEdges);
2376  m_corrHists["ChargeAsymmetry"].fillCorrHistsX(hit, hit.chargeAsymmetry);
2377  m_corrHists["ChargeLRplus"].fillCorrHistsX(hit, hit.chargeLRplus);
2378  m_corrHists["ChargeLRminus"].fillCorrHistsX(hit, hit.chargeLRminus);
2379  m_corrHists["SOverN"].fillCorrHistsX(hit, hit.sOverN);
2380  m_corrHists["WidthProj"].fillCorrHistsX(hit, hit.projWidth);
2381  m_corrHists["WidthDiff"].fillCorrHistsX(hit, hit.projWidth - static_cast<float>(hit.widthX));
2382 
2383  sector.WidthVsWidthProjected->Fill(hit.projWidth, hit.widthX);
2384  sector.PWidthVsWidthProjected->Fill(hit.projWidth, hit.widthX);
2385 
2386  sector.WidthDiffVsMaxStrip->Fill(hit.maxStrip, hit.projWidth - static_cast<float>(hit.widthX));
2387  sector.PWidthDiffVsMaxStrip->Fill(hit.maxStrip, hit.projWidth - static_cast<float>(hit.widthX));
2388 
2389  sector.WidthDiffVsSigmaXHit->Fill(hit.errXHit, hit.projWidth - static_cast<float>(hit.widthX));
2390  sector.PWidthDiffVsSigmaXHit->Fill(hit.errXHit, hit.projWidth - static_cast<float>(hit.widthX));
2391 
2392  sector.WidthVsPhiSensX->Fill(hit.phiSensX * 180. / M_PI, hit.widthX);
2393  sector.PWidthVsPhiSensX->Fill(hit.phiSensX * 180. / M_PI, hit.widthX);
2394  }
2395 
2396  // Hit Parameters
2397  m_corrHists["SigmaXHit"].fillCorrHistsX(hit, hit.errXHit * 10000.);
2398  m_corrHists["SigmaXTrk"].fillCorrHistsX(hit, hit.errXTrk * 10000.);
2399  m_corrHists["SigmaX"].fillCorrHistsX(hit, hit.errX * 10000.);
2400 
2401  m_corrHists["PhiSens"].fillCorrHistsX(hit, hit.phiSens * 180. / M_PI);
2402  m_corrHists["PhiSensX"].fillCorrHistsX(hit, hit.phiSensX * 180. / M_PI);
2403  m_corrHists["PhiSensY"].fillCorrHistsX(hit, hit.phiSensY * 180. / M_PI);
2404 
2405  sector.XHit->Fill(hit.xHit);
2406  sector.XTrk->Fill(hit.xTrk);
2407  sector.SigmaX2->Fill(hit.errX2 * 10000. * 10000.);
2408 
2409  sector.ResX->Fill(hit.resX * 10000.);
2410  sector.NorResX->Fill(hit.norResX);
2411 
2412  sector.ProbX->Fill(hit.probX);
2413 
2414  sector.PhiSensXVsBarycentreX->Fill(hit.baryStripX, hit.phiSensX * 180. / M_PI);
2415  sector.PPhiSensXVsBarycentreX->Fill(hit.baryStripX, hit.phiSensX * 180. / M_PI);
2416 }
2417 
2419  TrackerSectorStruct& sector) {
2420  std::map<std::string, TrackerSectorStruct::CorrelationHists>& m_corrHists(sector.m_correlationHistsY);
2421  // Do not fill anything for strip
2422  if (!hit.isPixelHit)
2423  return;
2424 
2425  // Cluster Parameters
2426  m_corrHists["WidthY"].fillCorrHistsY(hit, hit.widthY);
2427  m_corrHists["BaryStripY"].fillCorrHistsY(hit, hit.baryStripY);
2428 
2429  m_corrHists["ChargePixel"].fillCorrHistsY(hit, hit.chargePixel);
2430  m_corrHists["ClusterProbXY"].fillCorrHistsY(hit, hit.clusterProbabilityXY);
2431  m_corrHists["ClusterProbQ"].fillCorrHistsY(hit, hit.clusterProbabilityQ);
2432  m_corrHists["ClusterProbXYQ"].fillCorrHistsY(hit, hit.clusterProbabilityXYQ);
2433  m_corrHists["LogClusterProb"].fillCorrHistsY(hit, hit.logClusterProbability);
2434  m_corrHists["IsOnEdge"].fillCorrHistsY(hit, hit.isOnEdge);
2435  m_corrHists["HasBadPixels"].fillCorrHistsY(hit, hit.hasBadPixels);
2436  m_corrHists["SpansTwoRoc"].fillCorrHistsY(hit, hit.spansTwoRoc);
2437  m_corrHists["QBin"].fillCorrHistsY(hit, hit.qBin);
2438 
2439  // Hit Parameters
2440  m_corrHists["SigmaYHit"].fillCorrHistsY(hit, hit.errYHit * 10000.);
2441  m_corrHists["SigmaYTrk"].fillCorrHistsY(hit, hit.errYTrk * 10000.);
2442  m_corrHists["SigmaY"].fillCorrHistsY(hit, hit.errY * 10000.);
2443 
2444  m_corrHists["PhiSens"].fillCorrHistsY(hit, hit.phiSens * 180. / M_PI);
2445  m_corrHists["PhiSensX"].fillCorrHistsY(hit, hit.phiSensX * 180. / M_PI);
2446  m_corrHists["PhiSensY"].fillCorrHistsY(hit, hit.phiSensY * 180. / M_PI);
2447 
2448  sector.YHit->Fill(hit.yHit);
2449  sector.YTrk->Fill(hit.yTrk);
2450  sector.SigmaY2->Fill(hit.errY2 * 10000. * 10000.);
2451 
2452  sector.ResY->Fill(hit.resY * 10000.);
2453  sector.NorResY->Fill(hit.norResY);
2454 
2455  sector.ProbY->Fill(hit.probY);
2456 
2457  sector.PhiSensYVsBarycentreY->Fill(hit.baryStripY, hit.phiSensY * 180. / M_PI);
2458  sector.PPhiSensYVsBarycentreY->Fill(hit.baryStripY, hit.phiSensY * 180. / M_PI);
2459 }
2460 
2462  unsigned int goodHitsPerTrack(trackStruct.v_hitParams.size());
2463 
2464  if (parameterSet_.getParameter<bool>("applyTrackCuts")) {
2465  // which tracks to take? need min. nr. of selected hits?
2466  if (goodHitsPerTrack < minGoodHitsPerTrack_)
2467  return;
2468  }
2469 
2470  for (auto const& i_hit : trackStruct.v_hitParams) {
2471  // Put here from earlier method
2472  if (i_hit.hitState == TrackStruct::notAssignedToSectors)
2473  continue;
2474 
2475  for (auto& i_sector : m_tkSector_) {
2476  bool moduleInSector(false);
2477  for (auto const& i_hitSector : i_hit.v_sector) {
2478  if (i_sector.first == i_hitSector) {
2479  moduleInSector = true;
2480  break;
2481  }
2482  }
2483  if (!moduleInSector)
2484  continue;
2485 
2486  if (!calculateApe_)
2487  continue;
2488 
2489  if (i_hit.goodXMeasurement) {
2490  for (auto const& i_errBins : m_resErrBins_) {
2491  // Separate the bins for residual resolution w/o APE, to be consistent within iterations where APE will change (have same hit always in same bin)
2492  // So also fill this value in the histogram sigmaX
2493  // But of course use the normalized residual regarding the APE to have its influence in its width
2494  if (i_hit.errXWoApe < i_errBins.second.first || i_hit.errXWoApe >= i_errBins.second.second) {
2495  continue;
2496  }
2497  i_sector.second.m_binnedHists[i_errBins.first]["sigmaX"]->Fill(i_hit.errXWoApe);
2498  i_sector.second.m_binnedHists[i_errBins.first]["norResX"]->Fill(i_hit.norResX);
2499  break;
2500  }
2501  i_sector.second.ResX->Fill(i_hit.resX * 10000.);
2502  i_sector.second.NorResX->Fill(i_hit.norResX);
2503  }
2504 
2505  if (i_hit.goodYMeasurement) {
2506  for (auto const& i_errBins : m_resErrBins_) {
2507  // Separate the bins for residual resolution w/o APE, to be consistent within iterations where APE will change (have same hit always in same bin)
2508  // So also fill this value in the histogram sigmaY
2509  // But of course use the normalized residual regarding the APE to have its influence in its width
2510  if (i_hit.errYWoApe < i_errBins.second.first || i_hit.errYWoApe >= i_errBins.second.second) {
2511  continue;
2512  }
2513  i_sector.second.m_binnedHists[i_errBins.first]["sigmaY"]->Fill(i_hit.errYWoApe);
2514  i_sector.second.m_binnedHists[i_errBins.first]["norResY"]->Fill(i_hit.norResY);
2515  break;
2516  }
2517  i_sector.second.ResY->Fill(i_hit.resY * 10000.);
2518  i_sector.second.NorResY->Fill(i_hit.norResY);
2519  }
2520  }
2521  }
2522 }
2523 
2524 // -----------------------------------------------------------------------------------------------------------
2525 
2527  // Loop over sectors for calculating APE
2528  for (auto& i_sector : m_tkSector_) {
2529  // Loop over residual error bins to calculate APE for every bin
2530  for (auto const& i_errBins : i_sector.second.m_binnedHists) {
2531  std::map<std::string, TH1*> m_Hists = i_errBins.second;
2532 
2533  // Fitting Parameters
2534  double integralX = m_Hists["norResX"]->Integral();
2535  i_sector.second.EntriesX->SetBinContent(i_errBins.first, integralX);
2536 
2537  if (i_sector.second.isPixel) {
2538  double integralY = m_Hists["norResY"]->Integral();
2539  i_sector.second.EntriesY->SetBinContent(i_errBins.first, integralY);
2540  }
2541  }
2542  }
2543 }
2544 
2545 // -----------------------------------------------------------------------------------------------------------
2546 
2548  // we count SiStrip stereo modules as 2D if selected via countStereoHitAs2D_
2549  // (since they provide theta information)
2550  // --- NO, here it is always set to true ---
2551  if (!hit.isValid() || (hit.dimension() < 2 && !dynamic_cast<const SiStripRecHit1D*>(&hit))) {
2552  return false; // real RecHit1D - but SiStripRecHit1D depends on countStereoHitAs2D_
2553  } else {
2554  const DetId detId(hit.geographicalId());
2555  if (detId.det() == DetId::Tracker) {
2556  if (detId.subdetId() == PixelSubdetector::PixelBarrel || detId.subdetId() == PixelSubdetector::PixelEndcap) {
2557  return true; // pixel is always 2D
2558  } else { // should be SiStrip now
2559  const SiStripDetId stripId(detId);
2560  if (stripId.stereo())
2561  return true; // stereo modules
2562  else if (dynamic_cast<const SiStripRecHit1D*>(&hit) || dynamic_cast<const SiStripRecHit2D*>(&hit))
2563  return false; // rphi modules hit
2564  //the following two are not used any more since ages...
2565  else if (dynamic_cast<const SiStripMatchedRecHit2D*>(&hit))
2566  return true; // matched is 2D
2567  else if (dynamic_cast<const ProjectedSiStripRecHit2D*>(&hit)) {
2568  const ProjectedSiStripRecHit2D* pH = static_cast<const ProjectedSiStripRecHit2D*>(&hit);
2569  return (this->isHit2D(pH->originalHit())); // depends on original...
2570  } else {
2571  edm::LogError("UnkownType") << "@SUB=AlignmentTrackSelector::isHit2D"
2572  << "Tracker hit not in pixel, neither SiStripRecHit[12]D nor "
2573  << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
2574  return false;
2575  }
2576  }
2577  } else { // not tracker??
2578  edm::LogWarning("DetectorMismatch") << "@SUB=AlignmentTrackSelector::isHit2D"
2579  << "Hit not in tracker with 'official' dimension >=2.";
2580  return true; // dimension() >= 2 so accept that...
2581  }
2582  }
2583  // never reached...
2584 }
2585 
2586 // -----------------------------------------------------------------------------------------------------------
2587 
2588 // ------------ method called to for each event ------------
2591  edm::Handle<reco::BeamSpot> beamSpotHandle;
2592  iEvent.getByToken(offlinebeamSpot_, beamSpotHandle);
2593 
2594  if (beamSpotHandle.isValid()) {
2595  beamSpot = *beamSpotHandle;
2596  } else {
2597  edm::LogError("ApeEstimator") << "No beam spot available from EventSetup"
2598  << "\n...skip event";
2599  return;
2600  }
2601 
2603  iEvent.getByToken(tjTagToken_, m_TrajTracksMap);
2604 
2605  if (analyzerMode_)
2606  tkDetector_.TrkSize->Fill(m_TrajTracksMap->size());
2607 
2608  if (maxTracksPerEvent_ != 0 && m_TrajTracksMap->size() > maxTracksPerEvent_)
2609  return;
2610 
2611  //Creation of (traj,track)
2612  typedef std::pair<const Trajectory*, const reco::Track*> ConstTrajTrackPair;
2613  typedef std::vector<ConstTrajTrackPair> ConstTrajTrackPairCollection;
2614  ConstTrajTrackPairCollection trajTracks;
2615 
2617  for (i_trajTrack = m_TrajTracksMap->begin(); i_trajTrack != m_TrajTracksMap->end(); ++i_trajTrack) {
2618  trajTracks.push_back(ConstTrajTrackPair(&(*(*i_trajTrack).key), &(*(*i_trajTrack).val)));
2619  }
2620 
2621  //Loop over Tracks & Hits
2622  unsigned int trackSizeGood(0);
2623  ConstTrajTrackPairCollection::const_iterator iTrack;
2624  for (iTrack = trajTracks.begin(); iTrack != trajTracks.end(); ++iTrack) {
2625  const Trajectory* traj = (*iTrack).first;
2626  const reco::Track* track = (*iTrack).second;
2627 
2628  TrackStruct trackStruct;
2629  trackStruct.trkParams = this->fillTrackVariables(*track, *traj, beamSpot);
2630 
2631  if (trackCut_)
2632  continue;
2633 
2634  const std::vector<TrajectoryMeasurement> v_meas = (*traj).measurements();
2635 
2636  //Loop over Hits
2637  for (std::vector<TrajectoryMeasurement>::const_iterator i_meas = v_meas.begin(); i_meas != v_meas.end(); ++i_meas) {
2638  TrackStruct::HitParameterStruct hitParams = this->fillHitVariables(*i_meas, iSetup);
2639  if (this->hitSelected(hitParams))
2640  trackStruct.v_hitParams.push_back(hitParams);
2641  }
2642 
2643  if (analyzerMode_)
2644  this->fillHistsForAnalyzerMode(trackStruct);
2645  if (calculateApe_)
2646  this->fillHistsForApeCalculation(trackStruct);
2647 
2648  if (!trackStruct.v_hitParams.empty())
2649  ++trackSizeGood;
2650  }
2651  if (analyzerMode_ && trackSizeGood > 0)
2652  tkDetector_.TrkSizeGood->Fill(trackSizeGood);
2653 }
2654 
2655 // ------------ method called once each job just before starting event loop ------------
2657  this->hitSelection();
2658 
2659  this->sectorBuilder();
2660 
2661  this->residualErrorBinning();
2662 
2663  if (analyzerMode_)
2665 
2666  if (calculateApe_)
2668 
2669  if (analyzerMode_)
2670  this->bookTrackHists();
2671 }
2672 
2673 // ------------ method called once each job just after ending the event loop ------------
2675  if (calculateApe_)
2676  this->calculateAPE();
2677 
2678  edm::LogInfo("HitSelector") << "\nThere are " << counter1 << " negative Errors calculated\n";
2679 }
2680 
2681 //define this as a plug-in
ClusterRef cluster() const
virtual const Topology & topology() const =0
uint8_t maxCharge() const
double p() const
momentum vector magnitude
Definition: TrackBase.h:648
bool valid() const
Definition: LocalError.h:21
bool inDoubleInterval(const std::vector< double > &, const float) const
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
std::vector< unsigned int > v_sector
virtual float angularWidth() const =0
const edm::ParameterSet parameterSet_
float xx() const
Definition: LocalError.h:24
double d0Error() const
error on d0
Definition: TrackBase.h:847
virtual float length() const =0
std::map< std::string, std::vector< TH1 * > > m_sigmaX
void hitSelection()
float vv() const
unsigned short lost() const
Number of lost (=invalid) hits on track.
Definition: Track.h:204
float clusterProbability(unsigned int flags=0) const
T y() const
Definition: PV2DBase.h:46
ConstRecHitPointer const & recHit() const
void bookSectorHistsForAnalyzerMode()
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:630
uint32_t stereo() const
Definition: SiStripDetId.h:163
virtual const GeomDetType & type() const
Definition: GeomDet.cc:85
const_iterator end() const
last iterator over the map (read only)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::pair< uint16_t, uint16_t > chargeLR() const
virtual float originToIntersection() const =0
double theta() const
polar angle
Definition: TrackBase.h:612
std::map< std::string, std::vector< TH1 * > > m_sigmaY
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:124
LocalVector localDirection() const
ApeEstimator(const edm::ParameterSet &)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const bool calculateApe_
void sectorBuilder()
unsigned int counter1
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
float baryStrip() const
T y() const
Definition: PV3DBase.h:63
double etaError() const
error on eta
Definition: TrackBase.h:829
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:678
const Bounds & bounds() const
Definition: Surface.h:120
~ApeEstimator() override
int charge() const
PositionAndError2 radialPositionAndError2(const LocalPoint &, const LocalError &, const RadialStripTopology &)
unsigned int counter3
void fillHitHistsYForAnalyzerMode(const TrackStruct::HitParameterStruct &, TrackerSectorStruct &)
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
unsigned int counter4
uint16_t maxIndex() const
auto stripCharges() const -> decltype(cluster() ->amplitudes())
void calculateAPE()
bool isOnEdge() const
PositionAndError2 rectangularPositionAndError2(const LocalPoint &, const LocalError &)
TrackerDetectorStruct tkDetector_
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
LocalError positionError() const
const unsigned int maxTracksPerEvent_
bool hasBadPixels() const
TrackStruct::HitParameterStruct fillHitVariables(const TrajectoryMeasurement &, const edm::EventSetup &)
TProfile * PPhiSensYVsBarycentreY
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void residualErrorBinning()
const bool analyzerMode_
virtual float width() const =0
float localStripLength(const LocalPoint &) const override=0
DataContainer const & measurements() const
Definition: Trajectory.h:196
float signalOverNoise() const
float xy() const
Definition: LocalError.h:25
int iEvent
Definition: GenABIO.cc:224
float getLorentzAngle(const uint32_t &) const
TProfile * PPhiSensXVsBarycentreX
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:684
std::vector< HitParameterStruct > v_hitParams
uint16_t charge() const
float yy() const
Definition: LocalError.h:26
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:582
float stripAngle(float strip) const override=0
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:588
T sqrt(T t)
Definition: SSEVec.h:18
virtual int dimension() const =0
bool isHit2D(const TrackingRecHit &) const
TrackStruct::TrackParameterStruct fillTrackVariables(const reco::Track &, const Trajectory &, const reco::BeamSpot &)
double pt() const
track transverse momentum
Definition: TrackBase.h:654
TrackParameterStruct trkParams
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:808
edm::EDGetTokenT< reco::BeamSpot > offlinebeamSpot_
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
double phiError() const
error on phi
Definition: TrackBase.h:835
float uu() const
uint16_t width() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
bool hitSelected(TrackStruct::HitParameterStruct &) const
ClusterRef cluster() const
bool checkModulePositions(const float, const std::vector< double > &) const
void analyze(const edm::Event &, const edm::EventSetup &) override
T * make(const Args &...args) const
make new ROOT object
unsigned int counter6
const LocalTrajectoryError & localError() const
bool checkModuleBools(const bool, const std::vector< unsigned int > &) const
bool checkModuleIds(const unsigned int, const std::vector< unsigned int > &) const
virtual LocalPoint localPosition() const =0
std::map< unsigned int, TrackerSectorStruct > m_tkSector_
bool isValid() const
Definition: HandleBase.h:74
MeasurementError measurementError(const LocalPoint &, const LocalError &) const override=0
TProfile * PWidthVsWidthProjected
std::map< unsigned int, std::pair< double, double > > m_resErrBins_
bool checkIntervalsForSectors(const unsigned int sectorCounter, const std::vector< double > &) const
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
#define M_PI
void bookTrackHists()
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:642
double dzError() const
error on dz
Definition: TrackBase.h:859
SiStripRecHit2D originalHit() const
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
int qBin() const
Definition: DetId.h:18
virtual TrackingRecHit const * hit() const
virtual float detHeight() const =0
MeasurementPoint measurementPosition(const LocalPoint &) const override=0
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
uint16_t maxStrip() const
size_type size() const
map size
virtual float thickness() const =0
edm::Service< TFileService > fileService
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:479
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool checkModuleDirections(const int, const std::vector< int > &) const
Detector
Definition: DetId.h:26
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
bool isValid() const
bool inUintInterval(const std::vector< unsigned int > &, const unsigned int, const unsigned int=999) const
std::map< std::string, CorrelationHists > m_correlationHistsY
int sizeY() const
bool IsModuleUsable() const
std::map< unsigned int, ReducedTrackerTreeVariables > m_tkTreeVar_
Pixel cluster – collection of neighboring pixels above threshold.
void beginJob() override
fixed size matrix
unsigned int counter5
HLT enums.
virtual const GeomDetUnit * detUnit() const
void fillHistsForApeCalculation(const TrackStruct &)
unsigned short found() const
Number of valid hits on track.
Definition: Track.h:199
T get() const
Definition: EventSetup.h:71
unsigned int counter2
bool isPixel(HitType hitType)
bool spansTwoROCs() const
double y0() const
y coordinate
Definition: BeamSpot.h:66
const unsigned int minGoodHitsPerTrack_
float y() const
void setHitSelectionMap(const std::string &)
edm::EDGetTokenT< TrajTrackAssociationCollection > tjTagToken_
int charge() const
track electric charge
Definition: TrackBase.h:600
void fillHitHistsXForAnalyzerMode(const TrackStruct::HitParameterStruct &, TrackerSectorStruct &)
std::pair< const Trajectory *, const reco::Track * > ConstTrajTrackPair
Definition: Bounds.h:22
const_iterator begin() const
first iterator over the map (read only)
DetId geographicalId() const
void fillHistsForAnalyzerMode(const TrackStruct &)
void bookSectorHistsForApeCalculation()
virtual LocalError localPositionError() const =0
std::map< std::string, std::vector< double > > m_hitSelection_
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
void endJob() override
std::map< std::string, std::vector< unsigned int > > m_hitSelectionUInt_
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:624
T x() const
Definition: PV2DBase.h:45
TProfile * PWidthDiffVsSigmaXHit
int sizeX() const
std::pair< TrackStruct::HitState, PositionAndError2 > StatePositionAndError2
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
void setHitSelectionMapUInt(const std::string &)
T const * product() const
Definition: ESHandle.h:86
Definition: vlib.h:208
float x() const
LocalError const & localAlignmentError() const
Return local alligment error.
std::map< std::string, CorrelationHists > m_correlationHistsX
void statistics(const TrackerSectorStruct &, const int) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< unsigned int > v_rawId
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:45
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)
Our base class.
Definition: SiPixelRecHit.h:23
double x0() const
x coordinate
Definition: BeamSpot.h:64
StatePositionAndError2 positionAndError2(const LocalPoint &, const LocalError &, const TransientTrackingRecHit &)