CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FixTrackHitPattern.cc
Go to the documentation of this file.
3 
10 
11 // To convert detId to subdet/layer number.
14 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
21 
24 
28 
31 
33 
35 {
36  // Recalculate the inner and outer missing hit patterns. See header file for detailed comments.
37 
38  Result result;
39 
40  using namespace std;
41 
42  // Initialise Tracker geometry info (not sufficient to do this only on first call).
44  iSetup.get<TrackerRecoGeometryRecord>().get( tracker );
45 
46  // This is also needed to extrapolate amongst the tracker layers.
48  iSetup.get<NavigationSchoolRecord>().get("SimpleNavigationSchool",theSchool);
49  NavigationSetter junk(*theSchool);
50 
51  // This is needed to determine which sensors are functioning.
53  iSetup.get<CkfComponentsRecord>().get(measTk);
54  // Don't do this. It tries getting the tracker clusters, which do not exist in AOD.
55  // Hopefully not needed if one simply wants to know which sensors are active.
56  // measTk->update(iEvent);
57 
58  // Get the magnetic field and use it to define a propagator for extrapolating the track trajectory.
60  iSetup.get<IdealMagneticFieldRecord>().get(magField);
62 
63  // This is used to check if a track is compatible with crossing a sensor.
64  // Use +3.0 rather than default -3.0 here, so hit defined as inside acceptance if
65  // no more than 3*sigma outside detector edge, as opposed to more than 3*sigma inside detector edge.
66  Chi2MeasurementEstimator estimator(30.,3.0);
67 
68  // Get the track trajectory and detLayer at each valid hit on the track.
69  static GetTrackTrajInfo getTrackTrajInfo;
70  vector<GetTrackTrajInfo::Result> trkTrajInfo = getTrackTrajInfo.analyze(iSetup, track);
71  unsigned int nValidHits = trkTrajInfo.size();
72 
73  // Calculate the inner and outer hitPatterns on the track.
74  // i.e. The list of where the track should have had hits inside its innermost valid hit
75  // and outside its outmost valid hit.
76 
77  enum InnerOuter {INNER = 1, OUTER=2};
78  for (unsigned int inOut = INNER; inOut <= OUTER; inOut++) {
79 
80  // Get the track trajectory and detLayer at the inner/outermost valid hit.
81  unsigned int ihit = (inOut == INNER) ? 0 : nValidHits - 1;
82 
83  // Check that information about the track trajectory was available for this hit.
84  if (trkTrajInfo[ihit].valid) {
85  const DetLayer* detLayer = trkTrajInfo[ihit].detLayer;
86  const TrajectoryStateOnSurface& detTSOS = trkTrajInfo[ihit].detTSOS;
87  const FreeTrajectoryState* detFTS = detTSOS.freeTrajectoryState();
88 
89  // When getting inner hit pattern, must propagate track inwards from innermost layer with valid hit.
90  // When getting outer hit pattern, must propagate track outwards from outermost layer with valid hit.
91  const PropagationDirection direc = (inOut == INNER) ? oppositeToMomentum : alongMomentum;
92 
93  // Find all layers this track is compatible with in the desired direction, starting from this layer.
94  // Based on code in RecoTracker/TrackProducer/interface/TrackProducerBase.icc
95  // N.B. The following call uses code in RecoTracker/TkNavigation/src/SimpleBarrelNavigableLayer::compatibleLayers()
96  // and RecoTracker/TkNavigation/src/SimpleNavigableLayer::wellInside().
97  // These make some curious checks on the direction of the trajectory relative to its starting point,
98  // so care was required above when calculating detFTS.
99  vector<const DetLayer*> compLayers = detLayer->compatibleLayers(*detFTS, direc);
100  LogDebug("FTHP")<<"Number of inner/outer "<<inOut<<" layers intercepted by track = "<<compLayers.size()<<endl;
101 
102  int counter = 0;
103  reco::HitPattern newHitPattern;
104 
105  for(vector<const DetLayer *>::const_iterator it=compLayers.begin(); it!=compLayers.end();
106  ++it){
107  if ((*it)->basicComponents().empty()) {
108  //this should never happen. but better protect for it
109  edm::LogWarning("FixTrackHitPattern")<<"a detlayer with no components: I can not figure out a DetId from this layer. please investigate.";
110  continue;
111  }
112 
113  // Find the sensor in this layer which is closest to the track trajectory.
114  // And get the track trajectory at that sensor.
115  propagator.setPropagationDirection(direc);
116  vector< GeometricSearchDet::DetWithState > detWithState = (*it)->compatibleDets(detTSOS, propagator, estimator);
117  // Check that at least one sensor was compatible with the track trajectory.
118  if(detWithState.size() > 0) {
119  // Check that this sensor is functional
120  DetId id = detWithState.front().first->geographicalId();
121  const MeasurementDet* measDet = measTk->idToDet(id);
122  if(measDet->isActive()){
123  // Hence record that the track should have produced a hit here, but did not.
124  // Store the information in a HitPattern.
125  InvalidTrackingRecHit tmpHit(id, TrackingRecHit::missing);
126  newHitPattern.set(tmpHit, counter);
127  counter++;
128  } else {
129  // Missing hit expected here, since sensor was not functioning.
130  }
131  }
132  }//end loop over compatible layers
133 
134  // Store this result.
135  if (inOut == INNER) {
136  result.innerHitPattern = newHitPattern;
137  } else {
138  result.outerHitPattern = newHitPattern;
139  }
140 
141  // Print result for debugging.
142  LogDebug("FTHP")<<"Number of missing hits "<<newHitPattern.numberOfHits()<<"/"<<counter<<endl;
143  for (int j = 0; j < std::max(newHitPattern.numberOfHits(), counter); j++) {
144  uint32_t hp = newHitPattern.getHitPattern(j);
145  uint32_t subDet = newHitPattern.getSubStructure(hp);
146  uint32_t layer = newHitPattern.getLayer(hp);
147  uint32_t status = newHitPattern.getHitType(hp);
148  LogDebug("FTHP")<<" layer with no matched hit at counter="<<j<<" subdet="<<subDet<<" layer="<<layer<<" status="<<status<<endl;
149  }
150 
151  } else {
152  LogDebug("FTHP")<<"WARNING: could not calculate inner/outer hit pattern as trajectory info for inner/out hit missing"<<endl;
153  }
154 
155  }
156 
157  return result;
158 }
#define LogDebug(id)
static uint32_t getLayer(uint32_t pattern)
Definition: HitPattern.h:485
FreeTrajectoryState * freeTrajectoryState(bool withErrors=true) const
virtual bool isActive() const =0
static uint32_t getHitType(uint32_t pattern)
Definition: HitPattern.h:501
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
PropagationDirection
Result analyze(const edm::EventSetup &iSetup, const reco::Track &track)
const T & max(const T &a, const T &b)
tuple result
Definition: query.py:137
reco::HitPattern outerHitPattern
int numberOfHits() const
Definition: HitPattern.cc:213
static uint32_t getSubStructure(uint32_t pattern)
Definition: HitPattern.h:479
int j
Definition: DBlmapReader.cc:9
std::vector< Result > analyze(const edm::EventSetup &iSetup, const reco::Track &track)
Definition: DetId.h:20
reco::HitPattern innerHitPattern
const T & get() const
Definition: EventSetup.h:55
std::vector< const DetLayer * > compatibleLayers(NavigationDirection direction) const
Definition: DetLayer.cc:60
virtual void setPropagationDirection(PropagationDirection dir) const
Definition: Propagator.h:132
tuple status
Definition: ntuplemaker.py:245
uint32_t getHitPattern(int position) const
Definition: HitPattern.cc:144