CMS 3D CMS Logo

MuonCkfTrajectoryBuilder.cc
Go to the documentation of this file.
2 
4 
16 #include <sstream>
17 
19  : CkfTrajectoryBuilder(conf, iC),
20  theDeltaEta(conf.getParameter<double>("deltaEta")),
21  theDeltaPhi(conf.getParameter<double>("deltaPhi")),
22  theProximityPropagatorName(conf.getParameter<std::string>("propagatorProximity")),
23  theProximityPropagator(nullptr),
24  theEtaPhiEstimator(nullptr) {
25  //and something specific to me ?
26  theUseSeedLayer = conf.getParameter<bool>("useSeedLayer");
27  theRescaleErrorIfFail = conf.getParameter<double>("rescaleErrorIfFail");
28 }
29 
31 
33  CkfTrajectoryBuilder::setEvent_(iEvent, iSetup);
34 
35  edm::ESHandle<Propagator> propagatorProximityHandle;
36  iSetup.get<TrackingComponentsRecord>().get(theProximityPropagatorName, propagatorProximityHandle);
37  theProximityPropagator = propagatorProximityHandle.product();
38 
39  // theEstimator is set for this event in the base class
40  if (theEstimatorWatcher.check(iSetup) && theDeltaEta > 0 && theDeltaPhi > 0)
41  theEtaPhiEstimator = std::make_unique<EtaPhiEstimator>(theDeltaEta, theDeltaPhi, theEstimator);
42 }
43 
44 /*
45 std::string dumpMeasurement(const TrajectoryMeasurement & tm)
46 {
47  std::stringstream buffer;
48  buffer<<"layer pointer: "<<tm.layer()<<"\n"
49  <<"estimate: "<<tm.estimate()<<"\n"
50  <<"forward state: \n"
51  <<"x: "<<tm.forwardPredictedState().globalPosition()<<"\n"
52  <<"p: "<<tm.forwardPredictedState().globalMomentum()<<"\n"
53  <<"geomdet pointer from rechit: "<<tm.recHit()->det()<<"\n"
54  <<"detId: "<<tm.recHit()->geographicalId().rawId();
55  if (tm.recHit()->isValid()){
56  buffer<<"\n hit global x: "<<tm.recHit()->globalPosition()
57  <<"\n hit global error: "<<tm.recHit()->globalPositionError().matrix()
58  <<"\n hit local x:"<<tm.recHit()->localPosition()
59  <<"\n hit local error"<<tm.recHit()->localPositionError();
60  }
61  return buffer.str();
62 }
63 std::string dumpMeasurements(const std::vector<TrajectoryMeasurement> & v)
64 {
65  std::stringstream buffer;
66  buffer<<v.size()<<" total measurements\n";
67  for (std::vector<TrajectoryMeasurement>::const_iterator it = v.begin(); it!=v.end();++it){
68  buffer<<dumpMeasurement(*it);
69  buffer<<"\n";}
70  return buffer.str();
71 }
72 */
73 
75  const std::vector<const DetLayer*>& nl,
76  const TrajectoryStateOnSurface& currentState,
77  std::vector<TM>& result,
78  int& invalidHits,
79  const Propagator* prop) const {
80  for (std::vector<const DetLayer*>::const_iterator il = nl.begin(); il != nl.end(); il++) {
81  TSOS stateToUse = currentState;
82 
83  if (layer == (*il)) {
84  LogDebug("CkfPattern") << " self propagating in findCompatibleMeasurements.\n from: \n" << stateToUse;
85  //self navigation case
86  // go to a middle point first
88  GlobalPoint center(0, 0, 0);
89  stateToUse = middle.extrapolate(stateToUse, center, *prop);
90 
91  if (!stateToUse.isValid())
92  continue;
93  LogDebug("CkfPattern") << "to: " << stateToUse;
94  }
95 
97  std::vector<TM> tmp = layerMeasurements.measurements((**il), stateToUse, *prop, *theEstimator);
98 
99  if (tmp.size() == 1 && theEtaPhiEstimator) {
100  LogDebug("CkfPattern") << "only an invalid hit is found. trying differently";
101  tmp = layerMeasurements.measurements((**il), stateToUse, *prop, *theEtaPhiEstimator);
102  }
103  LogDebug("CkfPattern") << tmp.size() << " measurements returned by LayerMeasurements";
104 
105  if (!tmp.empty()) {
106  // FIXME durty-durty-durty cleaning: never do that please !
107  /* for (vector<TM>::iterator it = tmp.begin(); it!=tmp.end(); ++it)
108  {if (it->recHit()->det()==0) it=tmp.erase(it)--;}*/
109 
110  if (result.empty())
111  result = tmp;
112  else {
113  // keep one dummy TM at the end, skip the others
114  result.insert(result.end() - invalidHits, tmp.begin(), tmp.end());
115  }
116  invalidHits++;
117  }
118  }
119 
120  LogDebug("CkfPattern") << "starting from:\n"
121  << "x: " << currentState.globalPosition() << "\n"
122  << "p: " << currentState.globalMomentum() << "\n"
124 }
125 
127  const TempTrajectory& traj,
128  std::vector<TrajectoryMeasurement>& result) const {
129  int invalidHits = 0;
130 
131  std::vector<const DetLayer*> nl;
132 
133  if (traj.empty()) {
134  LogDebug("CkfPattern") << "using JR patch for no measurement case";
135  //what if there are no measurement on the Trajectory
136 
137  //set the currentState to be the one from the trajectory seed starting point
138  PTrajectoryStateOnDet ptod = seed.startingState();
139  DetId id(ptod.detId());
141  const Surface* surface = &g->surface();
142 
143  TrajectoryStateOnSurface currentState(
145 
146  //set the next layers to be that one the state is on
148 
149  if (theUseSeedLayer) {
150  {
151  //get the measurements on the layer first
152  LogDebug("CkfPattern") << "using the layer of the seed first.";
153  nl.push_back(l);
154  collectMeasurement(l, nl, currentState, result, invalidHits, theProximityPropagator);
155  }
156 
157  //if fails: try to rescale locally the state to find measurements
158  if ((unsigned int)invalidHits == result.size() && theRescaleErrorIfFail != 1.0 && !result.empty()) {
159  result.clear();
160  LogDebug("CkfPattern") << "using a rescale by " << theRescaleErrorIfFail << " to find measurements.";
161  TrajectoryStateOnSurface rescaledCurrentState = currentState;
162  rescaledCurrentState.rescaleError(theRescaleErrorIfFail);
163  invalidHits = 0;
164  collectMeasurement(l, nl, rescaledCurrentState, result, invalidHits, theProximityPropagator);
165  }
166  }
167 
168  //if fails: go to next layers
169  if (result.empty() || (unsigned int)invalidHits == result.size()) {
170  result.clear();
171  LogDebug("CkfPattern") << "Need to go to next layer to get measurements";
172  //the following will "JUMP" the first layer measurements
173  nl = theNavigationSchool->nextLayers(*l, *currentState.freeState(), traj.direction());
174  if (nl.empty()) {
175  LogDebug("CkfPattern") << " there was no next layer with wellInside. Use the next with no check.";
176  //means you did not get any compatible layer on the next 1/2 tracker layer.
177  // use the next layers with no checking
179  }
180  invalidHits = 0;
181  collectMeasurement(l, nl, currentState, result, invalidHits, forwardPropagator(seed));
182  }
183 
184  //if fails: this is on the next layers already, try rescaling locally the state
185  if (!result.empty() && (unsigned int)invalidHits == result.size() && theRescaleErrorIfFail != 1.0) {
186  result.clear();
187  LogDebug("CkfPattern") << "using a rescale by " << theRescaleErrorIfFail
188  << " to find measurements on next layers.";
189  TrajectoryStateOnSurface rescaledCurrentState = currentState;
190  rescaledCurrentState.rescaleError(theRescaleErrorIfFail);
191  invalidHits = 0;
192  collectMeasurement(l, nl, rescaledCurrentState, result, invalidHits, forwardPropagator(seed));
193  }
194 
195  } else //regular case
196  {
197  TSOS currentState(traj.lastMeasurement().updatedState());
198 
199  nl = theNavigationSchool->nextLayers(*traj.lastLayer(), *currentState.freeState(), traj.direction());
200  if (nl.empty()) {
201  LogDebug("CkfPattern") << " no next layers... going " << traj.direction() << "\n from: \n"
202  << currentState
203  << "\n from detId: " << traj.lastMeasurement().recHit()->geographicalId().rawId();
204  return;
205  }
206 
207  collectMeasurement(traj.lastLayer(), nl, currentState, result, invalidHits, forwardPropagator(seed));
208  }
209 
210  // sort the final result, keep dummy measurements at the end
211  if (result.size() > 1) {
212  sort(result.begin(), result.end() - invalidHits, TrajMeasLessEstim());
213  }
214 
215 #ifdef DEBUG_INVALID
216  bool afterInvalid = false;
217  for (std::vector<TM>::const_iterator i = result.begin(); i != result.end(); i++) {
218  if (!i->recHit().isValid())
219  afterInvalid = true;
220  if (afterInvalid && i->recHit().isValid()) {
221  edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!";
222  }
223  }
224 #endif
225 
226  //analyseMeasurements( result, traj);
227 }
#define LogDebug(id)
const Propagator * theProximityPropagator
T getParameter(std::string const &) const
const TrackingGeometry * geomTracker() const
ConstRecHitPointer const & recHit() const
bool empty() const
True if trajectory has no measurements.
const Propagator * forwardPropagator(const TrajectorySeed &seed) const
#define nullptr
GlobalPoint globalPosition() const
virtual const GeomDet * idToDet(DetId) const =0
static std::string dumpMeasurements(const std::vector< TrajectoryMeasurement > &v)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
const TrajectoryMeasurement & lastMeasurement() const
PropagationDirection direction() const
int iEvent
Definition: GenABIO.cc:224
void findCompatibleMeasurements(const TrajectorySeed &seed, const TempTrajectory &traj, std::vector< TrajectoryMeasurement > &result) const override
const std::string theProximityPropagatorName
FreeTrajectoryState const * freeState(bool withErrors=true) const
void setEvent_(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
const DetLayer * detLayer(const DetId &id) const
obsolete method. Use idToLayer() instead.
unsigned int detId() const
MuonCkfTrajectoryBuilder(const edm::ParameterSet &conf, edm::ConsumesCollector &iC)
const MeasurementTracker & measurementTracker() const
const MeasurementTrackerEvent * theMeasurementTracker
Definition: DetId.h:17
void collectMeasurement(const DetLayer *layer, const std::vector< const DetLayer * > &nl, const TrajectoryStateOnSurface &currentState, std::vector< TM > &result, int &invalidHits, const Propagator *) const
PTrajectoryStateOnDet const & startingState() const
std::unique_ptr< Chi2MeasurementEstimatorBase > theEtaPhiEstimator
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
std::vector< const DetLayer * > nextLayers(const DetLayer &detLayer, Args &&...args) const
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
edm::ESWatcher< BaseCkfTrajectoryBuilder::Chi2MeasurementEstimatorRecord > theEstimatorWatcher
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &fts, const GlobalPoint &vtx) const
extrapolation with default (=geometrical) propagator
GlobalVector globalMomentum() const
T get() const
Definition: EventSetup.h:73
const NavigationSchool * theNavigationSchool
TrajectoryStateOnSurface const & updatedState() const
void setEvent_(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
const DetLayer * lastLayer() const
Redundant method, returns the layer of lastMeasurement() .
tmp
align.sh
Definition: createJobs.py:716
T const * product() const
Definition: ESHandle.h:86
const GeometricSearchTracker * geometricSearchTracker() const
const Chi2MeasurementEstimatorBase * theEstimator