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