CMS 3D CMS Logo

DualTrajectoryFactory.cc
Go to the documentation of this file.
9 
10 #include <algorithm>
11 
14 
16 
18 public:
20  ~DualTrajectoryFactory() override;
21 
25  const reco::BeamSpot &beamSpot) const override;
26 
28  const ConstTrajTrackPairCollection &tracks,
30  const reco::BeamSpot &beamSpot) const override;
31 
32  DualTrajectoryFactory *clone() const override { return new DualTrajectoryFactory(*this); }
33 
34 protected:
39  };
40 
42 
44  const Surface &surface,
45  const MagneticField *magField) const;
46 
47  double theMass;
48 };
49 
53 
55  : TrajectoryFactoryBase(config), theMass(config.getParameter<double>("ParticleMass")) {
56  edm::LogInfo("Alignment") << "@SUB=DualTrajectoryFactory"
57  << "mass: " << theMass;
58 }
59 
61 
65 
67  setup.get<IdealMagneticFieldRecord>().get(magneticField);
68  if (magneticField->inTesla(GlobalPoint(0., 0., 0.)).mag2() < 1.e-6) {
69  edm::LogWarning("Alignment") << "@SUB=DualTrajectoryFactory::trajectories"
70  << "B-field in z is " << magneticField->inTesla(GlobalPoint(0., 0., 0.)).z()
71  << ": You should probably use the DualBzeroTrajectoryFactory\n"
72  << "or fix this code to switch automatically as the ReferenceTrajectoryFactory does.";
73  }
74 
75  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
76 
77  while (itTracks != tracks.end()) {
78  const DualTrajectoryInput input = this->referenceStateAndRecHits(*itTracks);
79  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
80  if (input.refTsos.isValid()) {
82  config.useBeamSpot = useBeamSpot_;
83  config.includeAPEs = includeAPEs_;
86  input.refTsos, input.fwdRecHits, input.bwdRecHits, magneticField.product(), beamSpot, config));
87  trajectories.push_back(ptr);
88  }
89 
90  ++itTracks;
91  }
92 
93  return trajectories;
94 }
95 
97  const edm::EventSetup &setup,
100  const reco::BeamSpot &beamSpot) const {
102 
103  if (tracks.size() != external.size()) {
104  edm::LogInfo("ReferenceTrajectories")
105  << "@SUB=DualTrajectoryFactory::trajectories"
106  << "Inconsistent input:\n"
107  << "\tnumber of tracks = " << tracks.size() << "\tnumber of external predictions = " << external.size();
108  return trajectories;
109  }
110 
112  setup.get<IdealMagneticFieldRecord>().get(magneticField);
113  if (magneticField->inTesla(GlobalPoint(0., 0., 0.)).mag2() < 1.e-6) {
114  edm::LogWarning("Alignment") << "@SUB=DualTrajectoryFactory::trajectories"
115  << "B-field in z is " << magneticField->inTesla(GlobalPoint(0., 0., 0.)).z()
116  << ": You should probably use the DualBzeroTrajectoryFactory\n"
117  << "or fix this code to switch automatically as the ReferenceTrajectoryFactory does.";
118  }
119 
120  ConstTrajTrackPairCollection::const_iterator itTracks = tracks.begin();
121  ExternalPredictionCollection::const_iterator itExternal = external.begin();
122 
123  while (itTracks != tracks.end()) {
125  // Check input: If all hits were rejected, the TSOS is initialized as invalid.
126  if (input.refTsos.isValid()) {
127  if ((*itExternal).isValid()) {
128  TrajectoryStateOnSurface propExternal =
129  propagateExternal(*itExternal, input.refTsos.surface(), magneticField.product());
130 
131  if (!propExternal.isValid())
132  continue;
133 
135  config.useBeamSpot = useBeamSpot_;
136  config.includeAPEs = includeAPEs_;
139  propExternal, input.fwdRecHits, input.bwdRecHits, magneticField.product(), beamSpot, config));
140 
141  AlgebraicSymMatrix externalParamErrors(asHepMatrix<5>(propExternal.localError().matrix()));
142  ptr->setParameterErrors(externalParamErrors);
143  trajectories.push_back(ptr);
144  } else {
146  config.useBeamSpot = useBeamSpot_;
147  config.includeAPEs = includeAPEs_;
150  input.refTsos, input.fwdRecHits, input.bwdRecHits, magneticField.product(), beamSpot, config));
151  trajectories.push_back(ptr);
152  }
153  }
154 
155  ++itTracks;
156  ++itExternal;
157  }
158 
159  return trajectories;
160 }
161 
163  const ConstTrajTrackPair &track) const {
165 
166  // get the trajectory measurements in the correct order, i.e. reverse if needed
167  Trajectory::DataContainer allTrajMeas = this->orderedTrajectoryMeasurements(*track.first);
168  Trajectory::DataContainer usedTrajMeas;
169  Trajectory::DataContainer::iterator itM;
170  // get all relevant trajectory measurements
171  for (itM = allTrajMeas.begin(); itM != allTrajMeas.end(); itM++) {
172  if (useRecHit((*itM).recHit()))
173  usedTrajMeas.push_back(*itM);
174  }
175 
176  unsigned int iMeas = 0;
177  unsigned int nMeas = usedTrajMeas.size();
178  unsigned int nRefStateMeas = nMeas / 2;
179  // get the valid RecHits
180  for (itM = usedTrajMeas.begin(); itM != usedTrajMeas.end(); itM++, iMeas++) {
181  TransientTrackingRecHit::ConstRecHitPointer aRecHit = (*itM).recHit();
182 
183  if (iMeas < nRefStateMeas) {
184  input.bwdRecHits.push_back(aRecHit);
185  } else if (iMeas > nRefStateMeas) {
186  input.fwdRecHits.push_back(aRecHit);
187  } else { // iMeas == nRefStateMeas
188  if ((*itM).updatedState().isValid()) {
189  input.refTsos = (*itM).updatedState();
190  input.bwdRecHits.push_back(aRecHit);
191  input.fwdRecHits.push_back(aRecHit);
192  } else {
193  // if the tsos of the middle hit is not valid, try the next one ...
194  nRefStateMeas++;
195  input.bwdRecHits.push_back(aRecHit);
196  }
197  }
198  }
199 
200  // bring input.fwdRecHits into correct order
201  std::reverse(input.bwdRecHits.begin(), input.bwdRecHits.end());
202 
203  return input;
204 }
205 
207  const Surface &surface,
208  const MagneticField *magField) const {
210  const std::pair<TrajectoryStateOnSurface, double> tsosWithPath = propagator.propagateWithPath(external, surface);
211  return tsosWithPath.first;
212 }
213 
T mag2() const
Definition: PV3DBase.h:66
const TrajectoryStateOnSurface propagateExternal(const TrajectoryStateOnSurface &external, const Surface &surface, const MagneticField *magField) const
const DualTrajectoryInput referenceStateAndRecHits(const ConstTrajTrackPair &track) const
MaterialEffects materialEffects(void) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
virtual const Trajectory::DataContainer orderedTrajectoryMeasurements(const Trajectory &trajectory) const
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
AlignmentAlgorithmBase::ConstTrajTrackPair ConstTrajTrackPair
Definition: config.py:1
bool useRecHit(const TransientTrackingRecHit::ConstRecHitPointer &hitPtr) const
static std::string const input
Definition: EdmProvDump.cc:48
config
Definition: looper.py:291
TransientTrackingRecHit::ConstRecHitContainer fwdRecHits
const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const override
Produce the reference trajectories.
const SurfaceType & surface() const
std::vector< TrajectoryMeasurement > DataContainer
Definition: Trajectory.h:44
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T z() const
Definition: PV3DBase.h:64
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
std::vector< ConstRecHitPointer > ConstRecHitContainer
DualTrajectoryFactory(const edm::ParameterSet &config)
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
TransientTrackingRecHit::ConstRecHitContainer bwdRecHits
DualTrajectoryFactory * clone() const override
T get() const
Definition: EventSetup.h:71
PropagationDirection propagationDirection(void) const
CLHEP::HepSymMatrix AlgebraicSymMatrix
std::vector< TrajectoryStateOnSurface > ExternalPredictionCollection
#define DEFINE_EDM_PLUGIN(factory, type, name)
A factory that produces instances of class ReferenceTrajectory from a given TrajTrackPairCollection.
T const * product() const
Definition: ESHandle.h:86
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &fts, const Plane &plane) const override
propagation to plane with path length