CMS 3D CMS Logo

TwoBodyDecayTrajectoryFactory.cc
Go to the documentation of this file.
8 
12 
14 
17 
20 
22 
29 public:
33 
38 
42  const reco::BeamSpot &beamSpot) const override;
43 
47  const reco::BeamSpot &beamSpot) const override;
48 
49  TwoBodyDecayTrajectoryFactory *clone() const override { return new TwoBodyDecayTrajectoryFactory(*this); }
50 
51 protected:
53  const TwoBodyDecay &tbd,
54  const MagneticField *magField,
55  const reco::BeamSpot &beamSpot,
56  bool setParameterErrors) const;
57 
59 
61 
65 
70 };
71 
75 
78  : TrajectoryFactoryBase(config, 2, iC),
79  m_MagFieldToken(iC.esConsumes()),
80  m_globTackingToken(iC.esConsumes()),
81  theFitter(config) {
82  const edm::ParameterSet ppc = config.getParameter<edm::ParameterSet>("ParticleProperties");
83  thePrimaryMass = ppc.getParameter<double>("PrimaryMass");
84  thePrimaryWidth = ppc.getParameter<double>("PrimaryWidth");
85  theSecondaryMass = ppc.getParameter<double>("SecondaryMass");
86 
87  theNSigmaCutValue = config.getParameter<double>("NSigmaCut");
88  theChi2CutValue = config.getParameter<double>("Chi2Cut");
89  theUseRefittedStateFlag = config.getParameter<bool>("UseRefittedState");
90  theConstructTsosWithErrorsFlag = config.getParameter<bool>("ConstructTsosWithErrors");
91 }
92 
94 
98 
100  const GlobalTrackingGeometry *trackingGeometry = &setup.getData(m_globTackingToken);
101  if (tracks.size() == 2) {
102  // produce transient tracks from persistent tracks
103  std::vector<reco::TransientTrack> transientTracks(2);
104 
105  transientTracks[0] = reco::TransientTrack(*tracks[0].second, magneticField);
106  transientTracks[0].setTrackingGeometry(trackingGeometry);
107 
108  transientTracks[1] = reco::TransientTrack(*tracks[1].second, magneticField);
109  transientTracks[1].setTrackingGeometry(trackingGeometry);
110 
111  // estimate the decay parameters
113  TwoBodyDecay tbd = theFitter.estimate(transientTracks, vm);
114 
115  if (!tbd.isValid() || (tbd.chi2() > theChi2CutValue)) {
117  return trajectories;
118  }
119 
120  return constructTrajectories(tracks, tbd, magneticField, beamSpot, false);
121  } else {
122  edm::LogInfo("ReferenceTrajectories") << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
123  << "Need 2 tracks, got " << tracks.size() << ".\n";
124  }
125 
126  return trajectories;
127 }
128 
130  const edm::EventSetup &setup,
133  const reco::BeamSpot &beamSpot) const {
135 
137 
138  const GlobalTrackingGeometry *trackingGeometry = &setup.getData(m_globTackingToken);
139 
140  if (tracks.size() == 2 && external.size() == 2) {
141  if (external[0].isValid() && external[1].isValid()) // Include external estimates
142  {
143  // produce transient tracks from persistent tracks
144  std::vector<reco::TransientTrack> transientTracks(2);
145 
146  transientTracks[0] = reco::TransientTrack(*tracks[0].second, magneticField);
147  transientTracks[0].setTrackingGeometry(trackingGeometry);
148 
149  transientTracks[1] = reco::TransientTrack(*tracks[1].second, magneticField);
150  transientTracks[1].setTrackingGeometry(trackingGeometry);
151 
152  // estimate the decay parameters. the transient tracks are not really associated to the
153  // the external tsos, but this is o.k., because the only information retrieved from them
154  // is the magnetic field.
156  TwoBodyDecay tbd = theFitter.estimate(transientTracks, external, vm);
157 
158  if (!tbd.isValid() || (tbd.chi2() > theChi2CutValue)) {
160  return trajectories;
161  }
162 
163  return constructTrajectories(tracks, tbd, magneticField, beamSpot, true);
164  } else {
165  // Return without external estimate
166  trajectories = this->trajectories(setup, tracks, beamSpot);
167  }
168  } else {
169  edm::LogInfo("ReferenceTrajectories") << "@SUB=TwoBodyDecayTrajectoryFactory::trajectories"
170  << "Need 2 tracks, got " << tracks.size() << ".\n";
171  }
172 
173  return trajectories;
174 }
175 
178  const TwoBodyDecay &tbd,
179  const MagneticField *magField,
180  const reco::BeamSpot &beamSpot,
181  bool setParameterErrors) const {
183 
184  // get innermost valid trajectory state and hits from the tracks
187 
188  if (!(input1.first.isValid() && input2.first.isValid()))
189  return trajectories;
190 
191  // produce TwoBodyDecayTrajectoryState (input for TwoBodyDecayTrajectory)
192  TsosContainer tsos(input1.first, input2.first);
193  ConstRecHitCollection recHits(input1.second, input2.second);
194  TwoBodyDecayTrajectoryState trajectoryState(tsos, tbd, theSecondaryMass, magField);
195 
196  if (!trajectoryState.isValid()) {
198  return trajectories;
199  }
200 
201  // always use the refitted trajectory state for matching
202  // FIXME FIXME CLONE
203  //TrackingRecHit::ConstRecHitPointer updatedRecHit1( recHits.first.front()->clone( tsos.first ) );
204  // TrackingRecHit::ConstRecHitPointer updatedRecHit2( recHits.second.front()->clone( tsos.second ) );
205 
206  TrackingRecHit::ConstRecHitPointer updatedRecHit1(recHits.first.front());
207  TrackingRecHit::ConstRecHitPointer updatedRecHit2(recHits.second.front());
208 
209  bool valid1 = match(trajectoryState.trajectoryStates(true).first, updatedRecHit1);
210 
211  bool valid2 = match(trajectoryState.trajectoryStates(true).second, updatedRecHit2);
212 
213  if (!valid1 || !valid2) {
215  return trajectories;
216  }
217 
219  config.useBeamSpot = useBeamSpot_;
220  config.includeAPEs = includeAPEs_;
221  config.allowZeroMaterial = allowZeroMaterial_;
222  config.useRefittedState = theUseRefittedStateFlag;
223  config.constructTsosWithErrors = theConstructTsosWithErrorsFlag;
224  // set the flag for reversing the RecHits to false, since they are already in the correct order.
225  config.hitsAreReverse = false;
226  TwoBodyDecayTrajectory *result = new TwoBodyDecayTrajectory(trajectoryState, recHits, magField, beamSpot, config);
227  if (setParameterErrors && tbd.hasError())
228  result->setParameterErrors(tbd.covariance());
230  return trajectories;
231 }
232 
235  LocalPoint lp1 = state.localPosition();
236  LocalPoint lp2 = recHit->localPosition();
237 
238  double deltaX = lp1.x() - lp2.x();
239  double deltaY = lp1.y() - lp2.y();
240 
241  LocalError le = recHit->localPositionError();
242 
243  double varX = le.xx();
244  double varY = le.yy();
245 
246  return ((fabs(deltaX) / sqrt(varX) < theNSigmaCutValue) && (fabs(deltaY) / sqrt(varY) < theNSigmaCutValue));
247 }
248 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > m_globTackingToken
TwoBodyDecayTrajectory::ConstRecHitCollection ConstRecHitCollection
TwoBodyDecayVirtualMeasurement VirtualMeasurement
dictionary config
Read in AllInOne config in JSON format.
Definition: DMR_cfg.py:21
bool hasError(void) const
Definition: TwoBodyDecay.h:42
Definition: config.py:1
MaterialEffects materialEffects(void) const
#define input2
Definition: AMPTWrapper.h:159
U second(std::pair< T, U > const &p)
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const ReferenceTrajectoryCollection trajectories(const edm::EventSetup &setup, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot &beamSpot) const override
Produce the trajectories.
const ReferenceTrajectoryCollection constructTrajectories(const ConstTrajTrackPairCollection &tracks, const TwoBodyDecay &tbd, const MagneticField *magField, const reco::BeamSpot &beamSpot, bool setParameterErrors) const
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:19
#define input1
Definition: AMPTWrapper.h:139
TwoBodyDecayTrajectoryFactory * clone() const override
std::pair< TrajectoryStateOnSurface, TransientTrackingRecHit::ConstRecHitContainer > TrajectoryInput
Log< level::Info, false > LogInfo
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > m_MagFieldToken
auto const & tracks
cannot be loose
virtual const TrajectoryInput innermostStateAndRecHits(const ConstTrajTrackPair &track) const
double chi2(void) const
Definition: TwoBodyDecay.h:44
bool match(const TrajectoryStateOnSurface &state, const TransientTrackingRecHit::ConstRecHitPointer &recHit) const
const AlgebraicSymMatrix & covariance(void) const
Definition: TwoBodyDecay.h:34
std::pair< TrajectoryStateOnSurface, TrajectoryStateOnSurface > TsosContainer
TwoBodyDecayTrajectoryFactory(const edm::ParameterSet &config, edm::ConsumesCollector &iC)
bool isValid(void) const
Definition: TwoBodyDecay.h:46
std::vector< TrajectoryStateOnSurface > ExternalPredictionCollection
const TsosContainer & trajectoryStates(bool useRefittedState=true) const
PropagationDirection propagationDirection(void) const
#define DEFINE_EDM_PLUGIN(factory, type, name)
ReferenceTrajectoryBase::ReferenceTrajectoryPtr ReferenceTrajectoryPtr
virtual const TwoBodyDecay estimate(const std::vector< reco::TransientTrack > &tracks, const TwoBodyDecayVirtualMeasurement &vm) const
std::vector< ReferenceTrajectoryPtr > ReferenceTrajectoryCollection
TwoBodyDecayTrajectoryState::TsosContainer TsosContainer
std::pair< ConstRecHitContainer, ConstRecHitContainer > ConstRecHitCollection
float xx() const
Definition: LocalError.h:22