CMS 3D CMS Logo

DAFTrackProducerAlgorithm.cc
Go to the documentation of this file.
26 
28  : conf_(conf), minHits_(conf.getParameter<int>("MinHits")) {}
29 
31  const MagneticField* theMF,
33  const MeasurementTrackerEvent* measTk,
34  const TrajectoryFitter* theFitter,
35  const TransientTrackingRecHitBuilder* builder,
36  const MultiRecHitCollector* measurementCollector,
38  const reco::BeamSpot& bs,
39  AlgoProductCollection& algoResults,
40  TrajAnnealingCollection& trajann,
41  bool TrajAnnSaving_,
42  AlgoProductCollection& algoResultsBeforeDAF,
43  AlgoProductCollection& algoResultsAfterDAF) const {
44  LogDebug("DAFTrackProducerAlgorithm") << "Size of map: " << TTmap.size() << "\n";
45 
46  int cont = 0;
47  int nTracksChanged = 0;
48 
49  for (TrajTrackAssociationCollection::const_iterator itTTmap = TTmap.begin(); itTTmap != TTmap.end(); itTTmap++) {
50  const edm::Ref<std::vector<Trajectory> > BeforeDAFTraj = itTTmap->key;
51  std::vector<TrajectoryMeasurement> BeforeDAFTrajMeas = BeforeDAFTraj->measurements();
52  const reco::TrackRef BeforeDAFTrack = itTTmap->val;
53 
54  float ndof = 0;
55  Trajectory CurrentTraj;
56 
57  if (BeforeDAFTraj->isValid()) {
58  LogDebug("DAFTrackProducerAlgorithm") << "The trajectory #" << cont + 1 << " is valid. \n";
59 
60  //getting the MultiRecHit collection and the trajectory with a first fit-smooth round
61  std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> hits =
62  collectHits(*BeforeDAFTraj, measurementCollector, &*measTk);
63 
64  //new initial fit
65  CurrentTraj = fit(hits, theFitter, *BeforeDAFTraj);
66 
67  //starting the annealing program
68  for (std::vector<double>::const_iterator ian = updator->getAnnealingProgram().begin();
69  ian != updator->getAnnealingProgram().end();
70  ian++) {
71  if (CurrentTraj.isValid()) {
72  LogDebug("DAFTrackProducerAlgorithm") << "Seed direction is " << CurrentTraj.seed().direction()
73  << ".Traj direction is " << CurrentTraj.direction() << std::endl;
74 
75  //updating MultiRecHits and fit-smooth again
76  std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> curiterationhits =
77  updateHits(CurrentTraj, updator, &*measTk, *ian);
78  if (curiterationhits.first.size() < 3) {
79  LogDebug("DAFTrackProducerAlgorithm")
80  << "Rejecting trajectory with " << curiterationhits.first.size() << " hits" << std::endl;
81  CurrentTraj = Trajectory();
82  break;
83  }
84 
85  CurrentTraj = fit(curiterationhits, theFitter, CurrentTraj);
86 
87  //saving trajectory for each annealing cycle ...
88  if (TrajAnnSaving_) {
89  TrajAnnealing temp(CurrentTraj, *ian);
90  trajann.push_back(temp);
91  }
92 
93  LogDebug("DAFTrackProducerAlgorithm") << "done annealing value " << (*ian);
94 
95  } else
96  break;
97  } //end annealing cycle
98 
99  int percOfHitsUnchangedAfterDAF =
100  (1. * checkHits(*BeforeDAFTraj, CurrentTraj) / (1. * BeforeDAFTrajMeas.size())) * 100.;
101  LogDebug("DAFTrackProducerAlgorithm")
102  << "Ended annealing program with " << percOfHitsUnchangedAfterDAF << " unchanged." << std::endl;
103 
104  //computing the ndof keeping into account the weights
105  ndof = calculateNdof(CurrentTraj);
106 
107  //checking if the trajectory has the minimum number of valid hits ( weight (>1e-6) )
108  //in order to remove tracks with too many outliers.
109  int goodHits = countingGoodHits(CurrentTraj);
110 
111  if (goodHits >= minHits_) {
112  bool ok = buildTrack(CurrentTraj, algoResults, ndof, bs, &BeforeDAFTrack);
113  // or filtered?
114  if (ok)
115  cont++;
116 
117  //saving tracks before and after DAF
118  if ((100. - percOfHitsUnchangedAfterDAF) > 0.) {
119  bool okBefore = buildTrack(*BeforeDAFTraj, algoResultsBeforeDAF, ndof, bs, &BeforeDAFTrack);
120  bool okAfter = buildTrack(CurrentTraj, algoResultsAfterDAF, ndof, bs, &BeforeDAFTrack);
121  if (okBefore && okAfter)
122  nTracksChanged++;
123  }
124  } else {
125  LogDebug("DAFTrackProducerAlgorithm") << "Rejecting trajectory with " << CurrentTraj.foundHits() << " hits";
126  }
127  } //end run on track collection
128  else {
129  LogDebug("DAFTrackProducerAlgorithm") << "Rejecting empty trajectory" << std::endl;
130  }
131  } //end run on track collection
132 
133  LogDebug("DAFTrackProducerAlgorithm") << "Number of Tracks found: " << cont << "\n";
134  LogDebug("DAFTrackProducerAlgorithm") << "Number of Tracks changed: " << nTracksChanged << "\n";
135 }
136 /*------------------------------------------------------------------------------------------------------*/
137 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> DAFTrackProducerAlgorithm::collectHits(
138  const Trajectory vtraj,
139  const MultiRecHitCollector* measurementCollector,
140  const MeasurementTrackerEvent* measTk) const {
141  LogDebug("DAFTrackProducerAlgorithm") << "Calling DAFTrackProducerAlgorithm::collectHits";
142 
143  //getting the traj measurements from the MeasurementCollector
145  std::vector<TrajectoryMeasurement> collectedmeas = measurementCollector->recHits(vtraj, measTk);
146 
147  //if the MeasurementCollector is empty, make an "empty" pair
148  //else taking the collected measured hits and building the pair
149  if (collectedmeas.empty())
151 
152  for (std::vector<TrajectoryMeasurement>::const_iterator iter = collectedmeas.begin(); iter != collectedmeas.end();
153  iter++) {
154  hits.push_back(iter->recHit());
155  }
156 
157  //TrajectoryStateWithArbitraryError() == Creates a TrajectoryState with the same parameters
158  // as the input one, but with "infinite" errors, i.e. errors so big that they don't
159  // bias a fit starting with this state.
160  //return std::make_pair(hits,TrajectoryStateWithArbitraryError()(collectedmeas.front().predictedState()));
161 
162  // I do not have to rescale the error because it is already rescaled in the fit code
163  TrajectoryStateOnSurface initialStateFromTrack = collectedmeas.front().predictedState();
164 
165  LogDebug("DAFTrackProducerAlgorithm")
166  << "Pair (hits, TSOS) with TSOS predicted(collectedmeas.front().predictedState())";
167  return std::make_pair(hits, initialStateFromTrack);
168 }
169 /*------------------------------------------------------------------------------------------------------*/
170 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> DAFTrackProducerAlgorithm::updateHits(
171  const Trajectory vtraj,
173  const MeasurementTrackerEvent* theMTE,
174  double annealing) const {
175  LogDebug("DAFTrackProducerAlgorithm") << "Calling DAFTrackProducerAlgorithm::updateHits";
177  std::vector<TrajectoryMeasurement> vmeas = vtraj.measurements();
178  std::vector<TrajectoryMeasurement>::reverse_iterator imeas;
179  unsigned int hitcounter = 1;
180 
181  //I run inversely on the trajectory obtained and update the state
182  for (imeas = vmeas.rbegin(); imeas != vmeas.rend(); imeas++, hitcounter++) {
183  DetId id = imeas->recHit()->geographicalId();
184  MeasurementDetWithData measDet = theMTE->idToDet(id);
185 
187  TrajectoryStateOnSurface combtsos;
188  if (hitcounter == vmeas.size())
189  combtsos = imeas->predictedState(); //fwd
190  else if (hitcounter == 1)
191  combtsos = imeas->backwardPredictedState(); //bwd
192  else
193  combtsos = combiner(imeas->backwardPredictedState(), imeas->predictedState());
194 
195  PrintHit(&*imeas->recHit(), combtsos);
196  if (imeas->recHit()->isValid()) {
197  TransientTrackingRecHit::RecHitPointer updated = updator->update(imeas->recHit(), combtsos, measDet, annealing);
198  hits.push_back(updated);
199  } else {
200  hits.push_back(imeas->recHit());
201  }
202  }
203 
204  TrajectoryStateOnSurface updatedStateFromTrack = vmeas.back().predictedState();
205 
206  //return std::make_pair(hits,TrajectoryStateWithArbitraryError()(vmeas.back().updatedState()));
207  LogDebug("DAFTrackProducerAlgorithm") << "Pair (hits, TSOS) with TSOS predicted (vmeas.back().predictedState())";
208 
209  return std::make_pair(hits, updatedStateFromTrack);
210 }
211 /*------------------------------------------------------------------------------------------------------*/
213  const std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>& hits,
214  const TrajectoryFitter* theFitter,
215  Trajectory vtraj) const {
216  //creating a new trajectory starting from the direction of the seed of the input one and the hits
217  Trajectory newVec = theFitter->fitOne(TrajectorySeed({}, {}, vtraj.seed().direction()), hits.first, hits.second);
218 
219  if (newVec.isValid())
220  return newVec;
221  else {
222  LogDebug("DAFTrackProducerAlgorithm") << "Fit no valid.";
223  return Trajectory();
224  }
225 }
226 /*------------------------------------------------------------------------------------------------------*/
228  AlgoProductCollection& algoResults,
229  float ndof,
230  const reco::BeamSpot& bs,
231  const reco::TrackRef* BeforeDAFTrack) const {
232  LogDebug("DAFTrackProducerAlgorithm") << " BUILDER " << std::endl;
233  ;
234  TrajectoryStateOnSurface innertsos;
235 
236  if (vtraj.isValid()) {
237  std::unique_ptr<Trajectory> theTraj(new Trajectory(vtraj));
238 
239  if (vtraj.direction() == alongMomentum) {
240  //if (theTraj->direction() == oppositeToMomentum) {
241  innertsos = vtraj.firstMeasurement().updatedState();
242  } else {
243  innertsos = vtraj.lastMeasurement().updatedState();
244  }
245 
246  TSCBLBuilderNoMaterial tscblBuilder;
247  TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(*(innertsos.freeState()), bs);
248 
249  if (tscbl.isValid() == false)
250  return false;
251 
253  math::XYZPoint pos(v.x(), v.y(), v.z());
255  math::XYZVector mom(p.x(), p.y(), p.z());
256 
257  // LogDebug("TrackProducer") <<v<<p<<std::endl;
258 
259  auto algo = (*BeforeDAFTrack)->algo();
260  std::unique_ptr<reco::Track> theTrack(new reco::Track(vtraj.chiSquared(),
261  ndof, //in the DAF the ndof is not-integer
262  pos,
263  mom,
264  tscbl.trackStateAtPCA().charge(),
266  algo));
267  theTrack->setQualityMask((*BeforeDAFTrack)->qualityMask());
268 
269  AlgoProduct aProduct{theTraj.release(), theTrack.release(), vtraj.direction(), 0};
270  algoResults.push_back(aProduct);
271 
272  return true;
273  } else {
274  LogDebug("DAFTrackProducerAlgorithm") << " BUILDER NOT POSSIBLE: traj is not valid" << std::endl;
275  ;
276  return false;
277  }
278 }
279 /*------------------------------------------------------------------------------------------------------*/
281  int ngoodhits = 0;
282  std::vector<TrajectoryMeasurement> vtm = traj.measurements();
283 
284  for (std::vector<TrajectoryMeasurement>::const_iterator tm = vtm.begin(); tm != vtm.end(); tm++) {
285  //if the rechit is valid
286  if (tm->recHit()->isValid()) {
287  SiTrackerMultiRecHit const& mHit = dynamic_cast<SiTrackerMultiRecHit const&>(*tm->recHit());
288  std::vector<const TrackingRecHit*> components = mHit.recHits();
289 
290  int iComp = 0;
291 
292  for (std::vector<const TrackingRecHit*>::const_iterator iter = components.begin(); iter != components.end();
293  iter++, iComp++) {
294  //if there is at least one component with weight higher than 1e-6 then the hit is not an outlier
295  if (mHit.weight(iComp) > 1e-6) {
296  ngoodhits++;
297  iComp++;
298  break;
299  }
300  }
301  }
302  }
303 
304  LogDebug("DAFTrackProducerAlgorithm") << "Original number of valid hits " << traj.foundHits()
305  << " -> hit with good weight (>1e-6) are " << ngoodhits;
306  return ngoodhits;
307 }
308 /*------------------------------------------------------------------------------------------------------*/
309 
311  std::vector<Trajectory>& input,
312  int minhits,
313  std::vector<Trajectory>& output,
314  const TransientTrackingRecHitBuilder* builder) const {
315  if (input.empty())
316  return;
317 
318  int ngoodhits = 0;
319  std::vector<TrajectoryMeasurement> vtm = input[0].measurements();
321 
322  //count the number of non-outlier and non-invalid hits
323  for (std::vector<TrajectoryMeasurement>::reverse_iterator tm = vtm.rbegin(); tm != vtm.rend(); tm++) {
324  //if the rechit is valid
325  if (tm->recHit()->isValid()) {
326  SiTrackerMultiRecHit const& mHit = dynamic_cast<SiTrackerMultiRecHit const&>(*tm->recHit());
327  std::vector<const TrackingRecHit*> components = mHit.recHits();
328  int iComp = 0;
329  bool isGood = false;
330  for (std::vector<const TrackingRecHit*>::const_iterator iter = components.begin(); iter != components.end();
331  iter++, iComp++) {
332  //if there is at least one component with weight higher than 1e-6 then the hit is not an outlier
333  if (mHit.weight(iComp) > 1e-6) {
334  ngoodhits++;
335  iComp++;
336  isGood = true;
337  break;
338  }
339  }
340  if (isGood) {
341  TkClonerImpl hc = static_cast<TkTransientTrackingRecHitBuilder const*>(builder)->cloner();
342  auto tempHit = hc.makeShared(tm->recHit(), tm->updatedState());
343  hits.push_back(tempHit);
344  } else
345  hits.push_back(std::make_shared<InvalidTrackingRecHit>(*tm->recHit()->det(), TrackingRecHit::missing));
346  } else {
347  TkClonerImpl hc = static_cast<TkTransientTrackingRecHitBuilder const*>(builder)->cloner();
348  auto tempHit = hc.makeShared(tm->recHit(), tm->updatedState());
349  hits.push_back(tempHit);
350  }
351  }
352 
353  LogDebug("DAFTrackProducerAlgorithm") << "Original number of valid hits " << input[0].foundHits()
354  << "; after filtering " << ngoodhits;
355  if (ngoodhits > input[0].foundHits())
356  edm::LogError("DAFTrackProducerAlgorithm")
357  << "Something wrong: the number of good hits from DAFTrackProducerAlgorithm::filter " << ngoodhits
358  << " is higher than the original one " << input[0].foundHits();
359 
360  if (ngoodhits < minhits)
361  return;
362 
363  TrajectoryStateOnSurface curstartingTSOS = input.front().lastMeasurement().updatedState();
364  LogDebug("DAFTrackProducerAlgorithm") << "starting tsos for final refitting " << curstartingTSOS;
365  //curstartingTSOS.rescaleError(100);
366 
367  output = fitter->fit(TrajectorySeed({}, {}, input.front().seed().direction()),
368  hits,
369  TrajectoryStateWithArbitraryError()(curstartingTSOS));
370 
371  LogDebug("DAFTrackProducerAlgorithm") << "After filtering " << output.size() << " trajectories";
372 }
373 /*------------------------------------------------------------------------------------------------------*/
375  if (!vtraj.isValid())
376  return 0;
377  float ndof = 0;
378  const std::vector<TrajectoryMeasurement>& meas = vtraj.measurements();
379  for (std::vector<TrajectoryMeasurement>::const_iterator iter = meas.begin(); iter != meas.end(); iter++) {
380  if (iter->recHit()->isValid()) {
381  SiTrackerMultiRecHit const& mHit = dynamic_cast<SiTrackerMultiRecHit const&>(*iter->recHit());
382  std::vector<const TrackingRecHit*> components = mHit.recHits();
383  int iComp = 0;
384  for (std::vector<const TrackingRecHit*>::const_iterator iter2 = components.begin(); iter2 != components.end();
385  iter2++, iComp++) {
386  if ((*iter2)->isValid())
387  ndof += ((*iter2)->dimension()) * mHit.weight(iComp);
388  }
389  }
390  }
391 
392  return ndof - 5;
393 }
394 //------------------------------------------------------------------------------------------------
395 int DAFTrackProducerAlgorithm::checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const {
396  std::vector<TrajectoryMeasurement> initmeasurements = iInitTraj.measurements();
397  std::vector<TrajectoryMeasurement> finalmeasurements = iFinalTraj.measurements();
398  std::vector<TrajectoryMeasurement>::iterator jmeas;
399  int nSame = 0;
400  int ihit = 0;
401 
402  if (initmeasurements.empty() || finalmeasurements.empty()) {
403  LogDebug("DAFTrackProducerAlgorithm") << "Initial or Final Trajectory empty.";
404  return 0;
405  }
406 
407  if (initmeasurements.size() != finalmeasurements.size()) {
408  LogDebug("DAFTrackProducerAlgorithm")
409  << "Initial Trajectory size(" << initmeasurements.size() << " hits) "
410  << "is different to final traj size (" << finalmeasurements.size() << ")! No checkHits possible! ";
411  return 0;
412  }
413 
414  for (jmeas = initmeasurements.begin(); jmeas != initmeasurements.end(); jmeas++) {
415  const TrackingRecHit* initHit = jmeas->recHit()->hit();
416  if (!initHit->isValid() && ihit == 0)
417  continue;
418 
419  if (initHit->isValid()) {
420  TrajectoryMeasurement imeas = finalmeasurements.at(ihit);
421  const TrackingRecHit* finalHit = imeas.recHit()->hit();
422  const TrackingRecHit* MaxWeightHit = nullptr;
423  float maxweight = 0;
424 
425  const SiTrackerMultiRecHit* mrh = dynamic_cast<const SiTrackerMultiRecHit*>(finalHit);
426  if (mrh) {
427  std::vector<const TrackingRecHit*> components = mrh->recHits();
428  std::vector<const TrackingRecHit*>::const_iterator icomp;
429  int hitcounter = 0;
430 
431  for (icomp = components.begin(); icomp != components.end(); icomp++) {
432  if ((*icomp)->isValid()) {
433  double weight = mrh->weight(hitcounter);
434  if (weight > maxweight) {
435  MaxWeightHit = *icomp;
436  maxweight = weight;
437  }
438  }
439  hitcounter++;
440  }
441  if (!MaxWeightHit)
442  continue;
443 
444  auto myref1 = reinterpret_cast<const BaseTrackerRecHit*>(initHit)->firstClusterRef();
445  auto myref2 = reinterpret_cast<const BaseTrackerRecHit*>(MaxWeightHit)->firstClusterRef();
446 
447  if (myref1 == myref2) {
448  nSame++;
449  } else {
450  LogDebug("DAFTrackProducerAlgorithm") << "diverso hit!" << std::endl;
451  TrajectoryStateOnSurface dummState;
452  LogTrace("DAFTrackProducerAlgorithm") << " This hit was:\n ";
453  PrintHit(initHit, dummState);
454  LogTrace("DAFTrackProducerAlgorithm") << " instead now is:\n ";
455  PrintHit(MaxWeightHit, dummState);
456  }
457  }
458  } else {
459  TrajectoryMeasurement imeas = finalmeasurements.at(ihit);
460  const TrackingRecHit* finalHit = imeas.recHit()->hit();
461  if (!finalHit->isValid()) {
462  nSame++;
463  }
464  }
465 
466  ihit++;
467  }
468 
469  return nSame;
470 }
471 
473  if (hit->isValid()) {
474  LogTrace("DAFTrackProducerAlgorithm")
475  << " Valid Hit with DetId " << hit->geographicalId().rawId() << " and dim:" << hit->dimension()
476  << " local position " << hit->localPosition() << " global position " << hit->globalPosition() << " and r "
477  << hit->globalPosition().perp();
478  if (tsos.isValid())
479  LogTrace("DAFTrackProducerAlgorithm") << " TSOS combtsos " << tsos.localPosition();
480 
481  } else {
482  LogTrace("DAFTrackProducerAlgorithm") << " Invalid Hit with DetId " << hit->geographicalId().rawId();
483  }
484 }
Trajectory fit(const std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > &hits, const TrajectoryFitter *theFitter, Trajectory vtraj) const
accomplishes the fitting-smoothing step for each annealing value
bool isValid() const
Definition: Trajectory.h:257
void PrintHit(const TrackingRecHit *const &hit, TrajectoryStateOnSurface &tsos) const
const CurvilinearTrajectoryError & curvilinearError() const
float calculateNdof(const Trajectory vtraj) const
bool buildTrack(const Trajectory, AlgoProductCollection &algoResults, float, const reco::BeamSpot &, const reco::TrackRef *) const
Construct Tracks to be put in the event.
std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > updateHits(const Trajectory vtraj, const SiTrackerMultiRecHitUpdator *updator, const MeasurementTrackerEvent *theMTE, double annealing) const
MeasurementDetWithData idToDet(const DetId &id) const
Previous MeasurementDetSystem interface.
std::vector< ConstRecHitPointer > RecHitContainer
std::pair< TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface > collectHits(const Trajectory vtraj, const MultiRecHitCollector *measurementCollector, const MeasurementTrackerEvent *measTk) const
virtual std::vector< TrajectoryMeasurement > recHits(const Trajectory &, const MeasurementTrackerEvent *theMTE) const =0
Definition: weight.py:1
float chiSquared() const
Definition: Trajectory.h:241
Log< level::Error, false > LogError
TrackingRecHit::ConstRecHitPointer makeShared(SiPixelRecHit const &hit, TrajectoryStateOnSurface const &tsos) const override
Definition: TkClonerImpl.cc:64
int foundHits() const
Definition: Trajectory.h:206
std::map< reco::TransientTrack, AlgebraicMatrix33 > TTmap
Definition: TTtoTTmap.h:11
key_type key() const
Accessor for product key.
Definition: Ref.h:250
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:150
#define LogTrace(id)
static std::string const input
Definition: EdmProvDump.cc:50
std::vector< TrajAnnealing > TrajAnnealingCollection
DataContainer const & measurements() const
Definition: Trajectory.h:178
GlobalPoint position() const
PropagationDirection direction() const
bool isValid() const
PropagationDirection const & direction() const
Definition: Trajectory.cc:133
TrackCharge charge() const
GlobalVector momentum() const
virtual Trajectory fitOne(const Trajectory &traj, fitType type=standard) const =0
void filter(const TrajectoryFitter *fitter, std::vector< Trajectory > &input, int minhits, std::vector< Trajectory > &output, const TransientTrackingRecHitBuilder *builder) const
std::shared_ptr< TrackingRecHit const > RecHitPointer
Definition: DetId.h:17
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:263
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
TrajectoryStateOnSurface const & updatedState() const
float weight(unsigned int i) const
int checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:166
int countingGoodHits(const Trajectory traj) const
typename Base::AlgoProductCollection AlgoProductCollection
FreeTrajectoryState const * freeState(bool withErrors=true) const
virtual TrackingRecHit const * hit() const
Definition: output.py:1
DAFTrackProducerAlgorithm(const edm::ParameterSet &conf)
void runWithCandidate(const TrackingGeometry *, const MagneticField *, const TrajTrackAssociationCollection &, const MeasurementTrackerEvent *measTk, const TrajectoryFitter *, const TransientTrackingRecHitBuilder *, const MultiRecHitCollector *measurementTracker, const SiTrackerMultiRecHitUpdator *, const reco::BeamSpot &, AlgoProductCollection &, TrajAnnealingCollection &, bool, AlgoProductCollection &, AlgoProductCollection &) const
Run the Final Fit taking TrackCandidates as input.
std::vector< Trajectory > fit(const Trajectory &traj, fitType type=standard) const
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
cont
load Luminosity info ##
Definition: generateEDF.py:620
ConstRecHitPointer const & recHit() const
#define LogDebug(id)